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 -------------------------------------------------------------------------------
12 License
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 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 defineTypeNameAndDebug(objectiveIncompressible, 0);
44 defineRunTimeSelectionTable(objectiveIncompressible, dictionary);
46 (
47  objective,
49  objective
50 );
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
54 objectiveIncompressible::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)
91 {
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  }
179  if (hasBoundarydJdTMVar1())
180  {
181  bdJdTMvar1Ptr_() *= oneOverNorm;
182  }
183  if (hasBoundarydJdTMVar2())
184  {
185  bdJdTMvar2Ptr_() *= oneOverNorm;
186  }
187  if (hasBoundarydJdnut())
188  {
189  bdJdnutPtr_() *= oneOverNorm;
190  }
191  if (hasBoundarydJdGradU())
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_() == dimensionedVector(dJdvPtr_().dimensions(), Zero);
244  }
245  if (hasdJdp())
246  {
247  dJdpPtr_() == dimensionedScalar(dJdpPtr_().dimensions(), Zero);
248  }
249  if (hasdJdT())
250  {
251  dJdTPtr_() == dimensionedScalar(dJdTPtr_().dimensions(), Zero);
252  }
253  if (hasdJdTMVar1())
254  {
255  dJdTMvar1Ptr_() ==
256  dimensionedScalar(dJdTMvar1Ptr_().dimensions(), Zero);
257  }
258  if (hasdJdTMVar2())
259  {
260  dJdTMvar2Ptr_() ==
261  dimensionedScalar(dJdTMvar2Ptr_().dimensions(), Zero);
262  }
263  if (hasBoundarydJdv())
264  {
265  bdJdvPtr_() == vector::zero;
266  }
267  if (hasBoundarydJdvn())
268  {
269  bdJdvnPtr_() == scalar(0);
270  }
271  if (hasBoundarydJdvt())
272  {
273  bdJdvtPtr_() == vector::zero;
274  }
275  if (hasBoundarydJdp())
276  {
277  bdJdpPtr_() == vector::zero;
278  }
279  if (hasBoundarydJdT())
280  {
281  bdJdTPtr_() == scalar(0);
282  }
283  if (hasBoundarydJdTMVar1())
284  {
285  bdJdTMvar1Ptr_() == scalar(0);
286  }
287  if (hasBoundarydJdTMVar2())
288  {
289  bdJdTMvar2Ptr_() == scalar(0);
290  }
291  if (hasBoundarydJdnut())
292  {
293  bdJdnutPtr_() == scalar(0);
294  }
295  if (hasBoundarydJdGradU())
296  {
297  bdJdGradUPtr_() == tensor::zero;
298  }
300  // Nullify geometric fields and sets nullified_ to true
302  }
303 }
304 
305 
307 {
308  // Figure out the availability of the adjoint turbulence model variables
309  // through their primal counterparts, since the contructor of the adjoint
310  // solver has not been completed yet
311  const incompressiblePrimalSolver& primSolver =
312  mesh_.lookupObject<incompressiblePrimalSolver>(primalSolverName_);
313  const autoPtr<incompressible::RASModelVariables>& rasVars =
314  primSolver.getIncoVars().RASModelVariables();
315 
316  if (rasVars().hasTMVar1())
317  {
318  const dimensionSet primalVarDims = rasVars->TMVar1Inst().dimensions();
319  const dimensionSet sourceDims = dimArea/pow3(dimTime)/primalVarDims;
320  dJdTMvar1Ptr_.reset
321  (
322  createZeroFieldPtr<scalar>
323  (
324  mesh_,
325  "ATMSource1",
326  sourceDims
327  )
328  );
329  }
330  if (rasVars().hasTMVar2())
331  {
332  const dimensionSet primalVarDims = rasVars->TMVar2Inst().dimensions();
333  const dimensionSet sourceDims = dimArea/pow3(dimTime)/primalVarDims;
334  dJdTMvar2Ptr_.reset
335  (
336  createZeroFieldPtr<scalar>
337  (
338  mesh_,
339  "ATMSource2",
340  sourceDims
341  )
342  );
343  }
344 }
345 
346 
348 (
349  const labelList& zoneIDs
350 ) const
351 {
352  label nCells(0);
353  for (const label zI : zoneIDs)
354  {
355  nCells += mesh_.cellZones()[zI].size();
356  }
357  reduce(nCells, sumOp<label>());
358  if (!nCells)
359  {
361  << "Provided cellZones include no cells"
362  << exit(FatalError);
363  }
364 }
365 
366 
368 (
369  autoPtr<volScalarField>& dJdTMvarPtr,
370  tmp<volScalarField>
371  (incompressibleAdjoint::adjointRASModel::*JacobianFunc)() const,
372  const volScalarField& JacobianMultiplier,
373  const labelList& zones
374 )
375 {
376  if (dJdTMvarPtr.good())
377  {
378  // nut Jacobians are currently computed in the adjoint turbulence
379  // models, though they would be better placed within the primal
380  // turbulence model.
381  // Safeguard the computation until the adjoint solver is complete
382  if (mesh_.foundObject<incompressibleAdjointSolver>(adjointSolverName_))
383  {
384  const incompressibleAdjointSolver& adjSolver =
385  mesh_.lookupObject<incompressibleAdjointSolver>
387  const autoPtr<incompressibleAdjoint::adjointRASModel>& adjointRAS =
388  adjSolver.getAdjointVars().adjointTurbulence();
389 
390  tmp<volScalarField> tnutJacobian = (adjointRAS->*JacobianFunc)();
391  const volScalarField& nutJacobian = tnutJacobian();
392 
393  volScalarField& dJdTMvar = dJdTMvarPtr();
394 
395  for (const label zI : zones)
396  {
397  const cellZone& zoneI = mesh_.cellZones()[zI];
398  for (const label cellI : zoneI)
399  {
400  dJdTMvar[cellI] =
401  JacobianMultiplier[cellI]*nutJacobian[cellI];
402  }
403  }
404  }
405  else
406  {
408  << "Skipping the computation of nutJacobian until "
409  << "the adjoint solver is complete"
410  << endl;
411  }
412  }
413 }
414 
415 
417 {
418  if (fieldNames_.found(matrix.psi().name()) && hasdJdv())
419  {
420  matrix += weight()*dJdv();
421  }
422 }
423 
424 
425 bool objectiveIncompressible::write(const bool valid) const
426 {
427  return objective::write(valid);
428 }
429 
430 
431 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
432 
433 } // End namespace Foam
434 
435 // ************************************************************************* //
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:53
virtual void update_meanValues()
Some objectives need to store some auxiliary values. If averaging is enabled, update these mean value...
dictionary dict
virtual void doNormalization()
Normalize all fields allocated by the objective.
static autoPtr< objectiveIncompressible > New(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
Return a reference to the selected turbulence model.
const labelIOList & zoneIDs
Definition: correctPhi.H:59
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...
const word adjointSolverName_
Definition: objective.H:65
const dictionary & dict() const
Return objective dictionary.
Definition: objective.C:91
autoPtr< volScalarField > dJdpPtr_
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
autoPtr< volScalarField > dJdTMvar2Ptr_
Second turbulence model variable.
autoPtr< boundaryScalarField > bdJdTMvar2Ptr_
Adjoint outlet turbulence model var 2.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
bool hasBoundarydJdTMVar1() const noexcept
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
virtual void update()
Update objective function derivatives.
autoPtr< volScalarField > dJdTPtr_
autoPtr< scalar > normFactor_
Normalization factor.
Definition: objective.H:101
void update_dJdTMvar(autoPtr< volScalarField > &dJdTMvarPtr, tmp< volScalarField >(incompressibleAdjoint::adjointRASModel::*JacobianFunc)() const, const volScalarField &JacobianMultiplier, const labelList &zones)
Compute dJdTMVar{1,2}.
virtual void doNormalization()
Normalize all fields allocated by the objective.
Definition: objective.C:359
bool hasdJdv() const noexcept
Inline functions for checking whether pointers are set or not.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
bool found(const T &val, label pos=0) const
Same as contains()
Definition: UList.H:879
void checkCellZonesSize(const labelList &zoneIDs) const
Check if cellZones provided include at least one cell.
autoPtr< boundaryTensorField > bdJdGradUPtr_
Term multiplying gradU variations.
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
Macros for easy insertion into run-time selection tables.
wordList fieldNames_
List of adjoint fields for which this objective will contribute sources to their equations.
Definition: objective.H:130
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
autoPtr< boundaryScalarField > bdJdnutPtr_
Jacobian wrt to nut.
virtual bool write(const bool valid=true) const
Write objective function history.
Definition: objective.C:523
dynamicFvMesh & mesh
bool computeMeanFields_
Definition: objective.H:68
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual void nullify()
Update objective function derivatives.
autoPtr< volVectorField > dJdvPtr_
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
const volVectorField & dJdv()
Contribution to field adjoint momentum eqs.
bool computed_
Whether the objective is computed or not.
Definition: objective.H:96
const fvMesh & mesh_
Definition: objective.H:63
autoPtr< boundaryVectorField > bdJdvPtr_
virtual bool write(const bool valid=true) const
Write objective function history.
virtual scalar weight() const
Return the objective function weight.
Definition: objective.C:347
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition: objective.H:55
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
autoPtr< boundaryVectorField > bdJdpPtr_
Adjoint (intlet,wall) velocity.
defineTypeNameAndDebug(combustionModel, 0)
bool hasBoundarydJdTMVar2() const noexcept
const word primalSolverName_
Definition: objective.H:66
autoPtr< boundaryScalarField > bdJdTPtr_
Adjoint outlet temperature.
void allocatedJdTurbulence()
Allocate fields related to the differentiation of turbulence models, if necessary.
void computeMeanFields()
Compute mean fields on the fly.
const incompressibleVars & vars_
virtual void addSource(fvVectorMatrix &matrix)
Vector sources can be given only to the adjoint momentum equations. Implemented in base objectiveInco...
dimensionedScalar pow3(const dimensionedScalar &ds)
virtual void nullify()
Nullify adjoint contributions.
Definition: objective.C:468
#define WarningInFunction
Report a warning using Foam::Warning.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual void update_dJdv()
Update vol and boundary fields and derivatives.
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:40
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:678
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
autoPtr< boundaryVectorField > bdJdvtPtr_
Adjoint outlet velocity.
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual void update()=0
Update objective function derivatives.
Definition: objective.C:447
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
autoPtr< boundaryScalarField > bdJdTMvar1Ptr_
Adjoint outlet turbulence model var 1.
Base class for primal incompressible solvers.
List< label > labelList
A List of labels.
Definition: List.H:62
autoPtr< boundaryScalarField > bdJdvnPtr_
Adjoint outlet pressure.
autoPtr< volScalarField > dJdTMvar1Ptr_
First turbulence model variable.
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:635
Abstract base class for objective functions in incompressible flows.
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127