Loading...
Searching...
No Matches
objectiveIncompressible.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
33#include "createZeroField.H"
35
36// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38namespace Foam
39{
40
41// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42
46(
47 objective,
48 objectiveIncompressible,
49 objective
50);
51
52// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53
54objectiveIncompressible::objectiveIncompressible
55(
56 const fvMesh& mesh,
57 const dictionary& dict,
58 const word& adjointSolverName,
59 const word& primalSolverName
60)
61:
62 objective(mesh, dict, adjointSolverName, primalSolverName),
63
64 vars_
65 (
66 mesh.lookupObject<incompressiblePrimalSolver>(primalSolverName).
67 getIncoVars()
68 ),
69
70 // Initialize pointers to nullptr.
71 // Not all of them are required for each objective function.
72 // Each child should allocate whatever is needed.
73
74 // Field adjoint Eqs
75 dJdvPtr_(nullptr),
76 dJdpPtr_(nullptr),
77 dJdTPtr_(nullptr),
78 dJdTMvar1Ptr_(nullptr),
79 dJdTMvar2Ptr_(nullptr),
80
81 // Adjoint boundary conditions
82 bdJdvPtr_(nullptr),
83 bdJdvnPtr_(nullptr),
84 bdJdvtPtr_(nullptr),
85 bdJdpPtr_(nullptr),
86 bdJdTPtr_(nullptr),
87 bdJdTMvar1Ptr_(nullptr),
88 bdJdTMvar2Ptr_(nullptr),
89 bdJdnutPtr_(nullptr),
90 bdJdGradUPtr_(nullptr)
92 computeMeanFields_ = vars_.computeMeanFields();
93}
94
95
96// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
97
99(
100 const fvMesh& mesh,
101 const dictionary& dict,
102 const word& adjointSolverName,
103 const word& primalSolverName
104)
105{
106 const word modelType(dict.get<word>("type"));
107
108 Info<< "Creating objective function : " << dict.dictName()
109 << " of type " << modelType << endl;
110
111 auto* ctorPtr = dictionaryConstructorTable(modelType);
112
113 if (!ctorPtr)
114 {
116 (
117 dict,
118 "objectiveIncompressible",
119 modelType,
120 *dictionaryConstructorTablePtr_
121 ) << exit(FatalIOError);
122 }
123
124 return autoPtr<objectiveIncompressible>
125 (
126 ctorPtr(mesh, dict, adjointSolverName, primalSolverName)
127 );
128}
129
130
131// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
132
134{
135 if (normalize_ && normFactor_)
136 {
137 const scalar oneOverNorm(1./normFactor_());
138
139 if (hasdJdv())
140 {
141 dJdvPtr_().primitiveFieldRef() *= oneOverNorm;
142 }
143 if (hasdJdp())
144 {
145 dJdpPtr_().primitiveFieldRef() *= oneOverNorm;
146 }
147 if (hasdJdT())
148 {
149 dJdTPtr_().primitiveFieldRef() *= oneOverNorm;
150 }
151 if (hasdJdTMVar1())
152 {
153 dJdTMvar1Ptr_().primitiveFieldRef() *= oneOverNorm;
154 }
155 if (hasdJdTMVar2())
156 {
157 dJdTMvar2Ptr_().primitiveFieldRef() *= oneOverNorm;
158 }
159 if (hasBoundarydJdv())
160 {
161 bdJdvPtr_() *= oneOverNorm;
162 }
163 if (hasBoundarydJdvn())
164 {
165 bdJdvnPtr_() *= oneOverNorm;
166 }
167 if (hasBoundarydJdvt())
168 {
169 bdJdvtPtr_() *= oneOverNorm;
170 }
171 if (hasBoundarydJdp())
172 {
173 bdJdpPtr_() *= oneOverNorm;
174 }
175 if (hasBoundarydJdT())
176 {
177 bdJdTPtr_() *= oneOverNorm;
178 }
180 {
181 bdJdTMvar1Ptr_() *= oneOverNorm;
182 }
184 {
185 bdJdTMvar2Ptr_() *= oneOverNorm;
186 }
187 if (hasBoundarydJdnut())
188 {
189 bdJdnutPtr_() *= oneOverNorm;
190 }
192 {
193 bdJdGradUPtr_() *= oneOverNorm;
194 }
196 // Normalize geometric fields
198 }
199}
200
201
203{
204 // Update geometric fields
206
207 // Update mean values here since they might be used in the
208 // subsequent functions
210
211 // volFields
212 update_dJdv();
213 update_dJdp();
214 update_dJdT();
217
218 // boundaryFields
228
229 // Divide everything with normalization factor
231
232 // Set objective as not computed, for the next optimisation cycle
233 computed_ = false;
234}
235
236
238{
239 if (!nullified_)
240 {
241 if (hasdJdv())
242 {
243 dJdvPtr_() == Zero;
244 }
245 if (hasdJdp())
246 {
247 dJdpPtr_() == Zero;
248 }
249 if (hasdJdT())
250 {
251 dJdTPtr_() == Zero;
252 }
253 if (hasdJdTMVar1())
254 {
255 dJdTMvar1Ptr_() == Zero;
256 }
257 if (hasdJdTMVar2())
258 {
259 dJdTMvar2Ptr_() == Zero;
260 }
261 if (hasBoundarydJdv())
262 {
263 bdJdvPtr_() == Zero;
264 }
265 if (hasBoundarydJdvn())
266 {
267 bdJdvnPtr_() == scalar(0);
268 }
269 if (hasBoundarydJdvt())
270 {
271 bdJdvtPtr_() == Zero;
272 }
273 if (hasBoundarydJdp())
274 {
275 bdJdpPtr_() == Zero;
276 }
277 if (hasBoundarydJdT())
278 {
279 bdJdTPtr_() == scalar(0);
280 }
282 {
283 bdJdTMvar1Ptr_() == scalar(0);
284 }
286 {
287 bdJdTMvar2Ptr_() == scalar(0);
288 }
289 if (hasBoundarydJdnut())
290 {
291 bdJdnutPtr_() == scalar(0);
292 }
294 {
295 bdJdGradUPtr_() == Zero;
296 }
298 // Nullify geometric fields and sets nullified_ to true
300 }
301}
302
303
305{
306 // Figure out the availability of the adjoint turbulence model variables
307 // through their primal counterparts, since the contructor of the adjoint
308 // solver has not been completed yet
309 const incompressiblePrimalSolver& primSolver =
310 mesh_.lookupObject<incompressiblePrimalSolver>(primalSolverName_);
311 const autoPtr<incompressible::RASModelVariables>& rasVars =
312 primSolver.getIncoVars().RASModelVariables();
313
314 if (rasVars().hasTMVar1())
315 {
316 const dimensionSet primalVarDims = rasVars->TMVar1Inst().dimensions();
317 const dimensionSet sourceDims = dimArea/pow3(dimTime)/primalVarDims;
318 dJdTMvar1Ptr_.reset
319 (
321 (
322 mesh_,
323 "ATMSource1",
324 sourceDims
325 )
326 );
327 }
328 if (rasVars().hasTMVar2())
329 {
330 const dimensionSet primalVarDims = rasVars->TMVar2Inst().dimensions();
331 const dimensionSet sourceDims = dimArea/pow3(dimTime)/primalVarDims;
332 dJdTMvar2Ptr_.reset
333 (
335 (
336 mesh_,
337 "ATMSource2",
338 sourceDims
339 )
340 );
341 }
342}
343
344
346(
347 const labelList& zoneIDs
348) const
349{
350 label nCells(0);
351 for (const label zI : zoneIDs)
352 {
353 nCells += mesh_.cellZones()[zI].size();
354 }
355 reduce(nCells, sumOp<label>());
356 if (!nCells)
357 {
359 << "Provided cellZones include no cells"
360 << exit(FatalError);
361 }
362}
363
364
366(
367 autoPtr<volScalarField>& dJdTMvarPtr,
368 tmp<volScalarField>
369 (incompressibleAdjoint::adjointRASModel::*JacobianFunc)() const,
370 const volScalarField& JacobianMultiplier,
371 const labelList& zones
372)
373{
374 if (dJdTMvarPtr.good())
375 {
376 // nut Jacobians are currently computed in the adjoint turbulence
377 // models, though they would be better placed within the primal
378 // turbulence model.
379 // Safeguard the computation until the adjoint solver is complete
380 if (mesh_.foundObject<incompressibleAdjointSolver>(adjointSolverName_))
381 {
382 const incompressibleAdjointSolver& adjSolver =
383 mesh_.lookupObject<incompressibleAdjointSolver>
385 const autoPtr<incompressibleAdjoint::adjointRASModel>& adjointRAS =
386 adjSolver.getAdjointVars().adjointTurbulence();
387
388 tmp<volScalarField> tnutJacobian = (adjointRAS->*JacobianFunc)();
389 const volScalarField& nutJacobian = tnutJacobian();
390
391 volScalarField& dJdTMvar = dJdTMvarPtr();
392
393 for (const label zI : zones)
394 {
395 const cellZone& zoneI = mesh_.cellZones()[zI];
396 for (const label cellI : zoneI)
397 {
398 dJdTMvar[cellI] =
399 JacobianMultiplier[cellI]*nutJacobian[cellI];
400 }
401 }
402 }
403 else
404 {
406 << "Skipping the computation of nutJacobian until "
407 << "the adjoint solver is complete"
408 << endl;
409 }
410 }
411}
412
413
415{
416 if (fieldNames_.found(matrix.psi().name()) && hasdJdv())
417 {
418 matrix += weight()*dJdv();
419 }
420}
421
422
423bool objectiveIncompressible::write(const bool valid) const
424{
425 return objective::write(valid);
426}
427
428
429// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
430
431} // End namespace Foam
432
433// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
bool good() const noexcept
True if the managed pointer is non-null.
Definition autoPtr.H:206
A subset of mesh cells.
Definition cellZone.H:61
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition fvMatrix.H:487
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Base class for incompressibleAdjoint solvers.
virtual const incompressibleAdjointVars & getAdjointVars() const
Access to the incompressible adjoint variables set.
const autoPtr< incompressibleAdjoint::adjointRASModel > & adjointTurbulence() const
Return const reference to the adjointRASModel.
Abstract base class for incompressible turbulence models.
Base class for primal incompressible solvers.
const incompressibleVars & getIncoVars() const
Access to the incompressible variables set.
const autoPtr< incompressible::RASModelVariables > & RASModelVariables() const
Return const reference to the turbulence model variables.
bool foundObject(const word &name, const bool recursive=false) const
Contains the named Type?
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
Abstract base class for objective functions in incompressible flows.
autoPtr< volScalarField > dJdTMvar2Ptr_
Second turbulence model variable.
static autoPtr< objectiveIncompressible > New(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
Return a reference to the selected turbulence model.
virtual bool write(const bool valid=true) const
Write objective function history.
bool hasdJdv() const noexcept
Inline functions for checking whether pointers are set or not.
autoPtr< boundaryTensorField > bdJdGradUPtr_
Term multiplying gradU variations.
autoPtr< boundaryScalarField > bdJdTMvar2Ptr_
Adjoint outlet turbulence model var 2.
autoPtr< volScalarField > dJdpPtr_
void update_dJdTMvar(autoPtr< volScalarField > &dJdTMvarPtr, tmp< volScalarField >(incompressibleAdjoint::adjointRASModel::*JacobianFunc)() const, const volScalarField &JacobianMultiplier, const labelList &zones)
Compute dJdTMVar{1,2}.
virtual void addSource(fvVectorMatrix &matrix)
Vector sources can be given only to the adjoint momentum equations. Implemented in base objectiveInco...
void allocatedJdTurbulence()
Allocate fields related to the differentiation of turbulence models, if necessary.
virtual void update_dJdv()
Update vol and boundary fields and derivatives.
virtual void nullify()
Update objective function derivatives.
const volVectorField & dJdv()
Contribution to field adjoint momentum eqs.
virtual void update_meanValues()
Some objectives need to store some auxiliary values. If averaging is enabled, update these mean value...
autoPtr< volScalarField > dJdTMvar1Ptr_
First turbulence model variable.
autoPtr< boundaryScalarField > bdJdTMvar1Ptr_
Adjoint outlet turbulence model var 1.
autoPtr< boundaryScalarField > bdJdnutPtr_
Jacobian wrt to nut.
autoPtr< boundaryScalarField > bdJdTPtr_
Adjoint outlet temperature.
autoPtr< boundaryScalarField > bdJdvnPtr_
Adjoint outlet pressure.
autoPtr< volScalarField > dJdTPtr_
autoPtr< boundaryVectorField > bdJdpPtr_
Adjoint (intlet,wall) velocity.
virtual void doNormalization()
Normalize all fields allocated by the objective.
virtual void update()
Update objective function derivatives.
autoPtr< boundaryVectorField > bdJdvtPtr_
Adjoint outlet velocity.
const incompressibleVars & vars_
autoPtr< boundaryVectorField > bdJdvPtr_
void checkCellZonesSize(const labelList &zoneIDs) const
Check if cellZones provided include at least one cell.
autoPtr< volVectorField > dJdvPtr_
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition objective.H:58
virtual bool write(const bool valid=true) const
Write objective function history.
Definition objective.C:527
const fvMesh & mesh_
Definition objective.H:63
virtual void update()=0
Update objective function derivatives.
Definition objective.C:453
wordList fieldNames_
List of adjoint fields for which this objective will contribute sources to their equations.
Definition objective.H:130
virtual scalar weight() const
Return the objective function weight.
Definition objective.C:353
virtual void nullify()
Nullify adjoint contributions.
Definition objective.C:474
autoPtr< scalar > normFactor_
Normalization factor.
Definition objective.H:101
bool computeMeanFields_
Definition objective.H:68
virtual void doNormalization()
Normalize all fields allocated by the objective.
Definition objective.C:365
const dictionary & dict() const
Return objective dictionary.
Definition objective.C:91
const word primalSolverName_
Definition objective.H:66
bool computed_
Whether the objective is computed or not.
Definition objective.H:96
const word adjointSolverName_
Definition objective.H:65
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition polyMesh.H:679
A class for managing temporary objects.
Definition tmp.H:75
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
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const labelIOList & zoneIDs
Definition correctPhi.H:59
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimArea(sqr(dimLength))
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
fvMatrix< vector > fvVectorMatrix
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
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
autoPtr< GeometricField< Type, fvPatchField, volMesh > > createZeroFieldPtr(const fvMesh &mesh, const word &name, const dimensionSet dims, bool printAllocation=false)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict