adjointRASModel.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 
30 #include "adjointRASModel.H"
31 #include "wallFvPatch.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace incompressibleAdjoint
39 {
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 defineTypeNameAndDebug(adjointRASModel, 0);
44 defineRunTimeSelectionTable(adjointRASModel, dictionary);
46 (
47  adjointTurbulenceModel,
50 );
51 
52 
53 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
54 
56 {
58  {
59  Info<< type() << "Coeffs" << coeffDict_ << endl;
60  }
61 }
62 
63 
65 {
66  const solverControl& solControl = adjointVars_.getSolverControl();
67  if (solControl.storeInitValues())
68  {
70  {
72  var1 == dimensionedScalar(var1.dimensions(), Zero);
73  }
74 
76  {
78  var2 == dimensionedScalar(var2.dimensions(), Zero);
79  }
80  }
81 }
82 
83 
85 {
86  const solverControl& solControl = adjointVars_.getSolverControl();
87  if (solControl.average())
88  {
90  {
92  (
93  new volScalarField
94  (
95  IOobject
96  (
97  getAdjointTMVariable1Inst().name() + "Mean",
98  mesh_.time().timeName(),
99  mesh_,
102  ),
104  )
105  );
106  }
107 
109  {
111  (
112  new volScalarField
113  (
114  IOobject
115  (
116  getAdjointTMVariable2Inst().name() + "Mean",
117  mesh_.time().timeName(),
118  mesh_,
121  ),
123  )
124  );
125  }
126  }
127 }
128 
129 
130 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
131 
132 adjointRASModel::adjointRASModel
133 (
134  const word& type,
135  incompressibleVars& primalVars,
137  objectiveManager& objManager,
138  const word& adjointTurbulenceModelName
139 )
140 :
141  adjointTurbulenceModel
142  (
143  primalVars,
144  adjointVars,
145  objManager,
146  adjointTurbulenceModelName
147  ),
149  (
150  IOobject
151  (
152  "adjointRASProperties",
153  primalVars.U().time().constant(),
154  primalVars.U().db(),
155  IOobject::MUST_READ_IF_MODIFIED,
156  IOobject::NO_WRITE
157  )
158  ),
159 
160  objectiveManager_(objManager),
161 
162  adjointTurbulence_(get<word>("adjointTurbulence")),
163  printCoeffs_(getOrDefault<Switch>("printCoeffs", false)),
164  coeffDict_(subOrEmptyDict(type + "Coeffs")),
165 
166  y_(mesh_),
167 
168  adjointTMVariable1Ptr_(nullptr),
169  adjointTMVariable2Ptr_(nullptr),
170  adjointTMVariablesBaseNames_(0),
171  adjointTMVariable1MeanPtr_(nullptr),
172  adjointTMVariable2MeanPtr_(nullptr),
173  adjMomentumBCSourcePtr_( createZeroBoundaryPtr<vector>(mesh_) ),
174  wallShapeSensitivitiesPtr_( createZeroBoundaryPtr<vector>(mesh_) ),
175  wallFloCoSensitivitiesPtr_( createZeroBoundaryPtr<vector>(mesh_) ),
176  includeDistance_(false),
177  changedPrimalSolution_(true)
178 {}
179 
180 
181 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
182 
184 (
185  incompressibleVars& primalVars,
187  objectiveManager& objManager,
188  const word& adjointTurbulenceModelName
189 )
190 {
191  const IOdictionary dict
192  (
193  IOobject
194  (
195  "adjointRASProperties",
196  primalVars.U().time().constant(),
197  primalVars.U().db(),
201  )
202  );
203 
204  const word modelType(dict.get<word>("adjointRASModel"));
205 
206  Info<< "Selecting adjointRAS turbulence model " << modelType << endl;
207 
208  auto* ctorPtr = dictionaryConstructorTable(modelType);
209 
210  if (!ctorPtr)
211  {
213  (
214  dict,
215  "adjointRASModel",
216  modelType,
217  *dictionaryConstructorTablePtr_
218  ) << exit(FatalIOError);
219  }
220 
222  (
223  ctorPtr
224  (
225  primalVars,
226  adjointVars,
227  objManager,
228  adjointTurbulenceModelName
229  )
230  );
231 }
232 
233 
234 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
235 
237 {
239 
241  {
242  y_.correct();
243  }
244 }
245 
246 
248 {
249  //if (regIOobject::read())
250 
251  // Bit of trickery : we are both IOdictionary ('adjointRASProperties') and
252  // an regIOobject from the adjointTurbulenceModel level. Problem is to
253  // distinguish between the two - we only want to reread the IOdictionary.
254 
255  bool ok = IOdictionary::readData
256  (
257  IOdictionary::readStream
258  (
260  )
261  );
263 
264  if (ok)
265  {
266  readEntry("adjointTurbulence", adjointTurbulence_);
267 
268  if (const dictionary* dictPtr = findDict(type() + "Coeffs"))
269  {
270  coeffDict_ <<= *dictPtr;
271  }
272 
273  return true;
274  }
275  else
276  {
277  return false;
278  }
279 }
280 
281 
283 {
285  {
286  // if pointer is not set, set it to a zero field
288  (
289  new volScalarField
290  (
291  IOobject
292  (
293  "adjointTMVariable1" + type(),
294  mesh_.time().timeName(),
295  mesh_,
298  ),
299  mesh_,
301  )
302  );
303  }
304 
305  return adjointTMVariable1Ptr_();
306 }
307 
308 
310 {
312  {
313  // if pointer is not set, set it to a zero field
315  (
316  new volScalarField
317  (
318  IOobject
319  (
320  "adjointTMVariable2" + type(),
321  mesh_.time().timeName(),
322  mesh_,
325  ),
326  mesh_,
328  )
329  );
330  }
331 
332  return adjointTMVariable2Ptr_();
333 }
334 
335 
337 {
339  {
341  }
342  else
343  {
344  return getAdjointTMVariable1Inst();
345  }
346 }
347 
348 
349 
351 {
353  {
355  }
356  else
357  {
358  return getAdjointTMVariable2Inst();
359  }
360 }
361 
364 {
365  return adjointTMVariable1Ptr_;
366 }
367 
370 {
371  return adjointTMVariable2Ptr_;
372 }
373 
376 {
378 }
379 
380 
382 {
383  dimensionSet dims
384  (
386  ? dimViscosity/adjointTMVariable1Ptr_().dimensions()
387  : dimless
388  );
389 
390  return
392  (
393  IOobject
394  (
395  "nutJacobianTMVar1"+type(),
396  mesh_.time().timeName(),
397  mesh_,
400  ),
401  mesh_,
402  dimensionedScalar(dims, Zero)
403  );
404 }
405 
406 
408 {
409  dimensionSet dims
410  (
412  ? dimViscosity/adjointTMVariable2Ptr_().dimensions()
413  : dimless
414  );
415 
416  return
418  (
419  IOobject
420  (
421  "nutJacobianTMVar2"+type(),
422  mesh_.time().timeName(),
423  mesh_,
426  ),
428  dimensionedScalar(dims, Zero)
429  );
430 }
431 
432 
434 (
435  tmp<volScalarField>& dNutdUMult
436 ) const
437 {
438  // Deliberately returning a null pointer in the base class
439  return nullptr;
440 }
441 
444 {
445  return tmp<scalarField>::New(mesh_.boundary()[patchI].size(), Zero);
446 }
447 
450 {
451  return tmp<scalarField>::New(mesh_.boundary()[patchI].size(), Zero);
452 }
453 
456 {
457  changedPrimalSolution_ = true;
458 }
459 
460 
462 {
463  const solverControl& solControl = adjointVars_.getSolverControl();
464  if (solControl.average())
465  {
467  {
470  }
472  {
475  }
476  }
477 }
478 
479 
481 {
482  const solverControl& solControl = adjointVars_.getSolverControl();
483  if (solControl.doAverageIter())
484  {
485  const label iAverageIter = solControl.averageIter();
486  scalar avIter(iAverageIter);
487  scalar oneOverItP1 = 1./(avIter+1);
488  scalar mult = avIter*oneOverItP1;
490  {
493  + getAdjointTMVariable1Inst()*oneOverItP1;
494  }
496  {
499  + getAdjointTMVariable2Inst()*oneOverItP1;
500  }
501  }
502 }
503 
504 
506 {
507  return includeDistance_;
508 }
509 
510 
511 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
512 
513 } // End namespace incompressible
514 } // End namespace Foam
515 
516 // ************************************************************************* //
dictionary coeffDict_
Model coefficients dictionary.
virtual tmp< scalarField > diffusionCoeffVar1(label patchI) const
Diffusion coefficient of the first primal and adjoint turbulence model equation. Needed for some adjo...
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
defineTypeNameAndDebug(adjointRASModel, 0)
dictionary dict
bool includeDistance_
Does the turbulence model include distances and should the adjoint to the distance field be computed...
type
Types of root.
Definition: Roots.H:52
autoPtr< volScalarField > adjointTMVariable1MeanPtr_
Adjoint turbulence model variable 1, mean value.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
nearWallDist y_
Near wall distance boundary field.
void computeMeanFields()
Average adjoint fields on the fly.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
defineRunTimeSelectionTable(adjointRASModel, dictionary)
Switch printCoeffs_
Flag to print the model coeffs at run-time.
static autoPtr< adjointRASModel > New(incompressibleVars &primalVars, incompressibleAdjointMeanFlowVars &adjointVars, objectiveManager &objManager, const word &adjointTurbulenceModelName=adjointTurbulenceModel::typeName)
Return a reference to the selected adjointRAS model.
Class for managing objective functions.
virtual bool read()
Read adjointRASProperties dictionary.
const dimensionSet dimViscosity
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
volScalarField & getAdjointTMVariable2()
Return non-constant reference to adjoint turbulence model variable 2.
virtual tmp< volVectorField > nutJacobianU(tmp< volScalarField > &dNutdUMult) const
Jacobian of nut wrt the flow velocity.
volScalarField & getAdjointTMVariable1()
Return non-constant reference to adjoint turbulence model variable 1.
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, any/none. Also accepts 0/1 as a string and shortcuts t/f, y/n.
Definition: Switch.H:77
void restoreInitValues()
Restore field values to the initial ones.
virtual void correct()=0
Solve the adjoint turbulence equations.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
Macros for easy insertion into run-time selection tables.
autoPtr< volScalarField > adjointTMVariable2MeanPtr_
Adjoint turbulence model variable 2, mean value.
autoPtr< volScalarField > & getAdjointTMVariable2InstPtr()
Return non-constant autoPtr to adjoint turbulence model variable 2.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
Switch adjointTurbulence_
Turbulence on/off flag.
void setChangedPrimalSolution()
Set flag of changed primal solution to true.
bool includeDistance() const
Should the adjoint to the eikonal equation be computed.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect...
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
Base class for solution control classes.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
bool changedPrimalSolution_
Has the primal solution changed?
virtual void correct()
Solve the adjoint turbulence equations.
autoPtr< typename GeometricField< Type, fvPatchField, volMesh >::Boundary > createZeroBoundaryPtr(const fvMesh &mesh, bool printAllocation=false)
void close()
Close Istream.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Abstract base class for incompressible turbulence models.
virtual tmp< scalarField > diffusionCoeffVar2(label patchI) const
Diffusion coefficient of the second primal and adjoint turbulence model equation. Needed for some adj...
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
volScalarField & getAdjointTMVariable1Inst()
Return non-constant reference to adjoint turbulence model variable 1.
Base class for solver control classes.
Definition: solverControl.H:45
Reading is optional [identical to LAZY_READ].
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:450
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
volScalarField & getAdjointTMVariable2Inst()
Return non-constant reference to adjoint turbulence model variable 2.
wordList adjointTMVariablesBaseNames_
Base names of the adjoint fields.
const word & name() const
Name function is needed to disambiguate those inherited from regIOobject and dictionary.
autoPtr< volScalarField > adjointTMVariable1Ptr_
Adjoint turbulence model variable 1.
autoPtr< volScalarField > adjointTMVariable2Ptr_
Adjoint turbulence model variable 2.
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition: IOobject.C:456
const volVectorField & U() const
Return const reference to velocity.
const wordList & getAdjointTMVariablesBaseNames() const
Return reference to the adjoint turbulence model variables base names.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
const word & constant() const noexcept
Return constant name.
Definition: TimePathsI.H:112
void resetMeanFields()
Reset mean fields to zero.
Abstract base class for incompressible adjoint turbulence models (RAS, LES and laminar).
const solverControl & getSolverControl() const
Return const reference to solverControl.
bool changing() const noexcept
Is mesh changing (topology changing and/or moving)
Definition: polyMesh.H:767
Manages the adjoint mean flow fields and their mean values.
addToRunTimeSelectionTable(adjointTurbulenceModel, adjointRASModel, adjointTurbulenceModel)
U
Definition: pEqn.H:72
bool average() const
Whether averaging is enabled or not.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Nothing to be read.
Automatically write from objectRegistry::writeObject()
messageStream Info
Information stream (stdout output on master, null elsewhere)
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
label & averageIter()
Return average iteration index reference.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
bool useAveragedFields() const
Use averaged fields? For solving the adjoint equations or computing sensitivities based on averaged f...
bool storeInitValues() const
Re-initialize.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:635
virtual tmp< volScalarField > nutJacobianTMVar2() const
Jacobian of nut wrt the second turbulence model variable.
virtual bool readData(Istream &)
The readData function required by regIOobject read operation.
virtual void correct()
Correct for mesh geom/topo changes.
Definition: nearWallDist.C:106
Do not request registration (bool: false)
autoPtr< volScalarField > & getAdjointTMVariable1InstPtr()
Return non-constant autoPtr to adjoint turbulence model variable 1.
Namespace for OpenFOAM.
const dimensionSet & dimensions() const noexcept
Return dimensions.
virtual tmp< volScalarField > nutJacobianTMVar1() const
Jacobian of nut wrt the first turbulence model variable.
virtual void printCoeffs()
Print model coefficients.
const dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary pointer if present (and a sub-dictionary) otherwise return nullptr...
Definition: dictionaryI.H:124
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