incompressibleAdjointSolver.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-2021 PCOpt/NTUA
9 Copyright (C) 2013-2021 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 "wallFvPatch.H"
34
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
42 (
46 );
47}
48
49
50// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51
53(
54 fvMesh& mesh,
55 const word& managerType,
56 const dictionary& dict,
57 const word& primalSolverName
58)
59:
60 adjointSolver(mesh, managerType, dict, primalSolverName),
61 primalVars_
62 (
63 mesh.lookupObjectRef<incompressiblePrimalSolver>(primalSolverName).
64 getIncoVars()
65 ),
66 ATCModel_(nullptr)
67{}
68
69
70// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
71
74(
75 fvMesh& mesh,
76 const word& managerType,
77 const dictionary& dict,
78 const word& primalSolverName
79)
80{
81 const word solverType(dict.get<word>("solver"));
82 auto* ctorPtr = dictionaryConstructorTable(solverType);
83
84 if (!ctorPtr)
85 {
87 (
88 dict,
89 "incompressibleAdjointSolver",
90 solverType,
91 *dictionaryConstructorTablePtr_
92 ) << exit(FatalIOError);
93 }
94
95 return
97 (
98 ctorPtr(mesh, managerType, dict, primalSolverName)
99 );
100}
101
102
103// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
104
106{
107 if (adjointSolver::readDict(dict))
108 {
109 return true;
110 }
111
112 return false;
113}
114
116{
117 return getAdjointVars().useSolverNameForFields();
118}
119
120
123{
124 return primalVars_;
125}
126
127
128
131{
132 const incompressibleAdjointVars& adjointVars =
133 refCast<incompressibleAdjointVars>(const_cast<variablesSet&>(vars_()));
134 return adjointVars;
135}
136
137
140{
141 incompressibleAdjointVars& adjointVars =
142 refCast<incompressibleAdjointVars>(const_cast<variablesSet&>(vars_()));
143 return adjointVars;
144}
145
146
147
150{
151 return ATCModel_;
152}
153
154
156{
157 return ATCModel_;
158}
159
160
162{
163 if (vars_)
164 {
165 getAdjointVars().adjointTurbulence()->setChangedPrimalSolution();
166 ATCModel_().updatePrimalBasedQuantities();
167 getAdjointVars().updatePrimalBasedQuantities();
168 }
169}
170
171
174{
175 /*
176 addProfiling
177 (
178 incompressibleAdjointSolver,
179 "incompressibleAdjointSolver::computeGradDxDbMultiplier"
180 );
181 */
183 (
184 getAdjointVars().adjointTurbulence()
185 );
186
187 const volScalarField& p = primalVars_.p();
188 const volVectorField& U = primalVars_.U();
189 const volScalarField& pa = getAdjointVars().pa();
190 const volVectorField& Ua = getAdjointVars().Ua();
191
192 // We only need to modify the boundaryField of gradU locally.
193 // If grad(U) is cached then
194 // a. The .ref() call fails since the tmp is initialised from a
195 // const ref
196 // b. we would be changing grad(U) for all other places in the code
197 // that need it
198 // So, always allocate new memory and avoid registering the new field
199 tmp<volTensorField> tgradU =
200 volTensorField::New("gradULocal", fvc::grad(U));
201 volTensorField& gradU = tgradU.ref();
203
204 // Explicitly correct the boundary gradient to get rid of
205 // the tangential component
206 forAll(mesh_.boundary(), patchI)
207 {
208 const fvPatch& patch = mesh_.boundary()[patchI];
209 if (isA<wallFvPatch>(patch))
210 {
211 tmp<vectorField> tnf = mesh_.boundary()[patchI].nf();
212 gradUbf[patchI] = tnf*U.boundaryField()[patchI].snGrad();
213 }
214 }
215
216 tmp<volScalarField> tnuEff = adjointRAS->nuEff();
217 tmp<volSymmTensorField> stress = tnuEff()*twoSymm(gradU);
218 // Note:
219 // term4 (Ua & grad(stress)) is numerically tricky. Its div leads to third
220 // order spatial derivs in E-SI based computations. Applying the product
221 // derivative rule (putting Ua inside the grad) gives better results in
222 // NACA0012, SA, WF. However, the original formulation should be kept at
223 // the boundary in order to respect the Ua boundary conditions (necessary
224 // for E-SI to give the same sens as FI). A mixed approach is hence
225 // followed
226
227 // Term 3, used also to allocated the return field
228 tmp<volTensorField> tgradUa = fvc::grad(Ua);
229 auto tflowTerm =
231 (
232 "flowTerm",
233 - tnuEff*(gradU & twoSymm(tgradUa()))
234 );
235 volTensorField& flowTerm = tflowTerm.ref();
236 // Term 4, only for the internal field
237 flowTerm.ref() +=
238 (
239 - (tgradUa & stress())
240 + fvc::grad(Ua & stress())
241 )().internalField();
242
243 // Boundary conditions from term 4
244 for (label idir = 0; idir < pTraits<vector>::nComponents; ++idir)
245 {
246 autoPtr<volVectorField> stressDirPtr
247 (
248 createZeroFieldPtr<vector>
249 (mesh_, "stressDir", stress().dimensions())
250 );
251 // Components need to be in the [0-5] range since stress is a
252 // volSymmTensorField
253 unzipRow(stress(), idir, stressDirPtr());
254 volTensorField gradStressDir(fvc::grad(stressDirPtr()));
255 forAll(mesh_.boundary(), pI)
256 {
257 if (!isA<coupledFvPatch>(mesh_.boundary()[pI]))
258 {
259 flowTerm.boundaryFieldRef()[pI] +=
260 Ua.component(idir)().boundaryField()[pI]
261 *gradStressDir.boundaryField()[pI];
262 }
263 }
264 }
265 // Release memory
266 stress.clear();
267
268 // Compute dxdb multiplier
269 flowTerm +=
270 // Term 1, ATC
271 ATCModel_->getFISensitivityTerm()
272 // Term 2
273 - fvc::grad(p)*Ua;
274
275 // Term 5
276 flowTerm += pa*tgradU;
277
278 // Term 6, from the adjoint turbulence model
279 flowTerm += T(adjointRAS->FISensitivityTerm());
280
281 // Term 7, term from objective functions
282 PtrList<objective>& functions
283 (objectiveManagerPtr_->getObjectiveFunctions());
284
285 for (objective& objI : functions)
286 {
287 if (objI.hasGradDxDbMult())
288 {
289 flowTerm += objI.weight()*objI.gradDxDbMultiplier();
290 }
291 }
292
293 flowTerm.correctBoundaryConditions();
294
295 //profiling::writeNow();
296
297 return (tflowTerm);
298}
299
300
302(
303 boundaryVectorField& sensitivityMap,
304 const labelHashSet& patchIDs,
305 const scalar dt
306)
307{
308 // Does nothing in base
309}
310
311
312// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Base class for adjoint solvers.
Definition: adjointSolver.H:60
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
Base class for incompressibleAdjoint solvers.
const incompressibleVars & getPrimalVars() const
Access to the incompressible primal variables set.
const autoPtr< ATCModel > & getATCModel() const
Access to the ATC model.
virtual void additionalSensitivityMapTerms(boundaryVectorField &sensitivityMap, const labelHashSet &patchIDs, const scalar dt)
virtual bool readDict(const dictionary &dict)
Read dict if updated.
virtual const incompressibleAdjointVars & getAdjointVars() const
Access to the incompressible adjoint variables set.
virtual bool useSolverNameForFields() const
Should solver name be appended to fields.
virtual tmp< volTensorField > computeGradDxDbMultiplier()
Compute the multiplier for grad(dxdb)
Class including all adjoint fields for incompressible flows.
Base class for primal incompressible solvers.
Base class for solution control classes.
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition: objective.H:62
A class for managing temporary objects.
Definition: tmp.H:65
void clear() const noexcept
Definition: tmpI.H:287
T & ref() const
Definition: tmpI.H:227
Base class for creating a set of variables.
Definition: variablesSet.H:50
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
volScalarField & p
const volScalarField & T
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Namespace for OpenFOAM.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
IOerror FatalIOError
void unzipRow(const FieldField< Field, SymmTensor< Cmpt > > &input, const direction idx, FieldField< Field, Vector< Cmpt > > &result)
Extract a symmTensor field field row (x,y,z) == (0,1,2)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333