Loading...
Searching...
No Matches
adjointSensitivity.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2007-2023 PCOpt/NTUA
9 Copyright (C) 2013-2023 FOSS GP
10 Copyright (C) 2019-2021 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
32#include "adjointSensitivity.H"
33#include "adjointSolver.H"
34#include "designVariables.H"
35#include "fvOptions.H"
36#include "reverseLinear.H"
37#include "sensitivity.H"
38
39// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41namespace Foam
42{
43
44// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46defineTypeNameAndDebug(adjointSensitivity, 0);
47defineRunTimeSelectionTable(adjointSensitivity, dictionary);
48
49
50// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51
52adjointSensitivity::adjointSensitivity
53(
54 const fvMesh& mesh,
55 const dictionary& dict,
57)
58:
59 sensitivity(mesh, dict),
61 derivatives_(0),
62 suffix_(word::null),
64 (
65 this->dict().getOrDefault<bool>
66 (
67 "includeDistance",
69 )
70 ),
71 eikonalSolver_(nullptr),
72 gradDxDbMult_(nullptr),
73 divDxDbMult_(nullptr),
74 dxdbMult_(nullptr),
75 dSfdbMult_(nullptr),
76 dnfdbMult_(nullptr),
77 dxdbDirectMult_(nullptr),
78 pointDxDbDirectMult_(nullptr),
79 bcDxDbMult_(nullptr),
80 optionsDxDbMult_(nullptr)
81{}
82
83
84// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
85
87(
88 const fvMesh& mesh,
89 const dictionary& dict,
91)
92{
93 const word sensType =
94 dict.optionalSubDict(mesh.name()).get<word>("sensitivityType");
95
96 Info<< "adjointSensitivity type : " << sensType << endl;
97
98 auto* ctorPtr = dictionaryConstructorTable(sensType);
99
100 if (!ctorPtr)
101 {
103 (
104 dict,
105 "adjointSensitivity",
106 sensType,
107 *dictionaryConstructorTablePtr_
108 ) << exit(FatalIOError);
109 }
110
111 return autoPtr<adjointSensitivity>
112 (
114 );
115}
116
117
118// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * //
119
121{
123 {
124 // The adjoint eikonal solver requires the parameterized patches
125 // as an argument, if they exist. Allocation will be managed by
126 // derived classes that have access to them
127 includeDistance_ = this->dict().getOrDefault<bool>
128 (
129 "includeDistance",
131 );
132
133 return true;
134 }
135
136 return false;
137}
138
141{
142 return false;
143}
144
145
147(
149)
150{
151 derivatives_ = designVars->assembleSensitivities(*this);
152}
153
154
156(
157 autoPtr<designVariables>& designVars
158)
160 assembleSensitivities(designVars);
161 write(type());
162 return derivatives_;
163}
164
167{
168 return derivatives_;
169}
170
171
173{
175 if (fieldSensPtr_)
176 {
177 fieldSensPtr_().primitiveFieldRef() = Zero;
178 }
179 if (eikonalSolver_)
180 {
181 eikonalSolver_->reset();
182 }
187}
188
189
190void adjointSensitivity::write(const word& baseName)
191{
192 sensitivity::write(baseName);
193}
194
195
196// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
197
198} // End namespace Foam
199
200// ************************************************************************* //
Abstract base class for adjoint-based sensitivities.
autoPtr< scalarField > divDxDbMult_
Multiplier of div(dx/db).
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
bool includeDistance_
Include distance variation in sensitivity computations.
autoPtr< boundaryVectorField > dxdbDirectMult_
Multiplier of dCf/db, found in the objective function.
virtual const scalarField & calculateSensitivities(autoPtr< designVariables > &designVars)
Calculates and returns sensitivity fields.
autoPtr< pointBoundaryVectorField > pointDxDbDirectMult_
Multiplier of dx/db computed at points, found in the objective function.
autoPtr< adjointEikonalSolver > eikonalSolver_
Adjoint eikonal equation solver.
autoPtr< boundaryVectorField > bcDxDbMult_
Multiplier of dx/db, coming from boundary conditions that depend on the geometry, like rotatingWallVe...
adjointSolver & adjointSolver_
Reference to the underlaying adjoint solver.
autoPtr< vectorField > optionsDxDbMult_
dx/db multiplier coming from fvOptions
virtual void write(const word &baseName=word::null)
Write sensitivity fields.
virtual bool readDict(const dictionary &dict)
Read dictionary if changed.
bool includeDistance() const
Should the adjoint eikonal PDE should be solved.
autoPtr< boundaryVectorField > dxdbMult_
Multiplier of face dx/db.
virtual bool computeDxDbInternalField() const
Should the parameterization compute the internalField of dxdb.
scalarField derivatives_
The sensitivity derivative values.
const scalarField & getSensitivities() const
Returns the sensitivity fields.
static autoPtr< adjointSensitivity > New(const fvMesh &mesh, const dictionary &dict, adjointSolver &adjointSolver)
Return a reference to the selected turbulence model.
autoPtr< boundaryVectorField > dnfdbMult_
Multiplier of dnf/db.
autoPtr< adjointMeshMovementSolver > adjointMeshMovementSolver_
Adjoint grid displacement solver.
autoPtr< boundaryVectorField > dSfdbMult_
Multiplier of dSf/db.
autoPtr< volTensorField > gradDxDbMult_
Multiplier of grad(dx/b).
virtual void assembleSensitivities(autoPtr< designVariables > &designVars)
Assemble sensitivities.
word suffix_
Append this word to files related to the sensitivities.
Base class for adjoint solvers.
virtual bool includeDistance() const
Does the adjoint to an equation computing distances need to taken into consideration.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
virtual void write(const word &baseName=word::null)
Write sensitivity fields.
Definition sensitivity.C:59
virtual bool readDict(const dictionary &dict)
Read dictionary if changed.
Definition sensitivity.C:51
const fvMesh & mesh() const
Return reference to mesh.
const dictionary & dict() const
Return the construction dictionary.
autoPtr< volScalarField > fieldSensPtr_
Definition sensitivity.H:74
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition error.H:637
Namespace for OpenFOAM.
messageStream Info
Information stream (stdout output on master, null elsewhere).
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
runTime write()
Macros to ease declaration of run-time selection tables.
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict