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-2019 PCOpt/NTUA
9  Copyright (C) 2013-2019 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.average())
68  {
70  {
72  (
73  new volScalarField
74  (
75  IOobject
76  (
77  getAdjointTMVariable1Inst().name() + "Mean",
78  mesh_.time().timeName(),
79  mesh_,
82  ),
84  )
85  );
86  }
87 
89  {
91  (
92  new volScalarField
93  (
94  IOobject
95  (
96  getAdjointTMVariable2Inst().name() + "Mean",
97  mesh_.time().timeName(),
98  mesh_,
101  ),
103  )
104  );
105  }
106  }
107 }
108 
109 
110 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
111 
112 adjointRASModel::adjointRASModel
113 (
114  const word& type,
115  incompressibleVars& primalVars,
117  objectiveManager& objManager,
118  const word& adjointTurbulenceModelName
119 )
120 :
121  adjointTurbulenceModel
122  (
123  primalVars,
124  adjointVars,
125  objManager,
126  adjointTurbulenceModelName
127  ),
129  (
130  IOobject
131  (
132  "adjointRASProperties",
133  primalVars.U().time().constant(),
134  primalVars.U().db(),
135  IOobject::MUST_READ_IF_MODIFIED,
136  IOobject::NO_WRITE
137  )
138  ),
139 
140  objectiveManager_(objManager),
141 
142  adjointTurbulence_(get<word>("adjointTurbulence")),
143  printCoeffs_(getOrDefault<Switch>("printCoeffs", false)),
144  coeffDict_(subOrEmptyDict(type + "Coeffs")),
145 
146  y_(mesh_),
147 
148  adjointTMVariable1Ptr_(nullptr),
149  adjointTMVariable2Ptr_(nullptr),
150  adjointTMVariablesBaseNames_(0),
151  adjointTMVariable1MeanPtr_(nullptr),
152  adjointTMVariable2MeanPtr_(nullptr),
153  adjMomentumBCSourcePtr_( createZeroBoundaryPtr<vector>(mesh_) ),
154  wallShapeSensitivitiesPtr_( createZeroBoundaryPtr<vector>(mesh_) ),
155  wallFloCoSensitivitiesPtr_( createZeroBoundaryPtr<vector>(mesh_) ),
156  includeDistance_(false),
157  changedPrimalSolution_(true)
158 {}
159 
160 
161 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
162 
164 (
165  incompressibleVars& primalVars,
167  objectiveManager& objManager,
168  const word& adjointTurbulenceModelName
169 )
170 {
171  const IOdictionary dict
172  (
173  IOobject
174  (
175  "adjointRASProperties",
176  primalVars.U().time().constant(),
177  primalVars.U().db(),
180  false // Do not register
181  )
182  );
183 
184  const word modelType(dict.get<word>("adjointRASModel"));
185 
186  Info<< "Selecting adjointRAS turbulence model " << modelType << endl;
187 
188  auto* ctorPtr = dictionaryConstructorTable(modelType);
189 
190  if (!ctorPtr)
191  {
193  (
194  dict,
195  "adjointRASModel",
196  modelType,
197  *dictionaryConstructorTablePtr_
198  ) << exit(FatalIOError);
199  }
200 
202  (
203  ctorPtr
204  (
205  primalVars,
206  adjointVars,
207  objManager,
208  adjointTurbulenceModelName
209  )
210  );
211 }
212 
213 
214 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
215 
217 {
219 
221  {
222  y_.correct();
223  }
224 }
225 
226 
228 {
229  //if (regIOobject::read())
230 
231  // Bit of trickery : we are both IOdictionary ('adjointRASProperties') and
232  // an regIOobject from the adjointTurbulenceModel level. Problem is to
233  // distinguish between the two - we only want to reread the IOdictionary.
234 
235  bool ok = IOdictionary::readData
236  (
237  IOdictionary::readStream
238  (
240  )
241  );
243 
244  if (ok)
245  {
246  readEntry("adjointTurbulence", adjointTurbulence_);
247 
248  if (const dictionary* dictPtr = findDict(type() + "Coeffs"))
249  {
250  coeffDict_ <<= *dictPtr;
251  }
252 
253  return true;
254  }
255  else
256  {
257  return false;
258  }
259 }
260 
261 
263 {
265  {
266  // if pointer is not set, set it to a zero field
268  (
269  new volScalarField
270  (
271  IOobject
272  (
273  "adjointTMVariable1" + type(),
274  mesh_.time().timeName(),
275  mesh_,
278  ),
279  mesh_,
281  )
282  );
283  }
284 
285  return adjointTMVariable1Ptr_();
286 }
287 
288 
290 {
292  {
293  // if pointer is not set, set it to a zero field
295  (
296  new volScalarField
297  (
298  IOobject
299  (
300  "adjointTMVariable2" + type(),
301  mesh_.time().timeName(),
302  mesh_,
305  ),
306  mesh_,
308  )
309  );
310  }
311 
312  return adjointTMVariable2Ptr_();
313 }
314 
315 
317 {
319  {
321  }
322  else
323  {
324  return getAdjointTMVariable1Inst();
325  }
326 }
327 
328 
329 
331 {
333  {
335  }
336  else
337  {
338  return getAdjointTMVariable2Inst();
339  }
340 }
341 
344 {
345  return adjointTMVariable1Ptr_;
346 }
347 
350 {
351  return adjointTMVariable2Ptr_;
352 }
353 
356 {
358 }
359 
360 
362 {
363  return
365  (
366  IOobject
367  (
368  "nutJacobianTMVar1"+type(),
369  mesh_.time().timeName(),
370  mesh_,
373  ),
374  mesh_,
376  (
378  Zero
379  )
380  );
381 }
382 
383 
385 {
386  return
388  (
389  IOobject
390  (
391  "nutJacobianTMVar2"+type(),
392  mesh_.time().timeName(),
393  mesh_,
396  ),
397  mesh_,
399  (
400  nut().dimensions()/adjointTMVariable2Ptr_().dimensions(),
402  )
403  );
404 }
405 
406 
408 (
409  tmp<volScalarField>& dNutdUMult
410 ) const
411 {
412  // Deliberately returning a null pointer in the base class
413  return nullptr;
414 }
415 
418 {
419  return tmp<scalarField>::New(mesh_.boundary()[patchI].size(), Zero);
420 }
421 
424 {
425  return tmp<scalarField>::New(mesh_.boundary()[patchI].size(), Zero);
426 }
427 
430 {
431  changedPrimalSolution_ = true;
432 }
433 
434 
436 {
437  const solverControl& solControl = adjointVars_.getSolverControl();
438  if (solControl.average())
439  {
441  {
444  }
446  {
449  }
450  }
451 }
452 
453 
455 {
456  const solverControl& solControl = adjointVars_.getSolverControl();
457  if (solControl.doAverageIter())
458  {
459  const label iAverageIter = solControl.averageIter();
460  scalar avIter(iAverageIter);
461  scalar oneOverItP1 = 1./(avIter+1);
462  scalar mult = avIter*oneOverItP1;
464  {
467  + getAdjointTMVariable1Inst()*oneOverItP1;
468  }
470  {
473  + getAdjointTMVariable2Inst()*oneOverItP1;
474  }
475  }
476 }
477 
478 
480 {
481  return includeDistance_;
482 }
483 
484 
485 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
486 
487 } // End namespace incompressible
488 } // End namespace Foam
489 
490 // ************************************************************************* //
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:120
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 incompressible objective functions.
virtual bool read()
Read adjointRASProperties dictionary.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
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
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:361
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:84
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:752
Base class for solution control classes.
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.
Reading required, file watched for runTime modification.
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 elements in the list.
Definition: UPtrListI.H:99
volScalarField & getAdjointTMVariable1Inst()
Return non-constant reference to adjoint turbulence model variable 1.
Base class for solver control classes.
Definition: solverControl.H:45
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:441
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:203
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 wordList & getAdjointTMVariablesBaseNames()
Return reference to the adjoint turbulence model variables base names.
const volVectorField & U() const
Return const reference to velocity.
virtual const volScalarField & nut() const
Return the turbulence viscosity.
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 time name of given scalar time formatted with the given precision.
Definition: Time.C:760
const word & constant() const noexcept
Return constant name.
Definition: TimePathsI.H:89
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:732
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.
const Time & time() const
Return Time associated with the objectRegistry.
Definition: IOobject.C:447
IOobject(const IOobject &)=default
Copy construct.
Nothing to be read.
Automatically write from objectRegistry::writeObject()
messageStream Info
Information stream (stdout output on master, null elsewhere)
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...
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:166
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:615
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
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:120
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:705
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:157