Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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
30#include "adjointRASModel.H"
31#include "wallFvPatch.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
39{
40
41// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42
46(
48 adjointRASModel,
49 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 }
73
75 {
77 }
78 }
79}
80
81
83{
84 const solverControl& solControl = adjointVars_.getSolverControl();
85 if (solControl.average())
86 {
88 {
90 (
92 (
93 IOobject
94 (
97 mesh_,
100 ),
102 )
103 );
104 }
105
107 {
109 (
111 (
112 IOobject
113 (
114 getAdjointTMVariable2Inst().name() + "Mean",
115 mesh_.time().timeName(),
116 mesh_,
119 ),
121 )
122 );
124 }
125}
126
127
128// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
129
130adjointRASModel::adjointRASModel
131(
132 const word& type,
133 incompressibleVars& primalVars,
134 incompressibleAdjointMeanFlowVars& adjointVars,
135 objectiveManager& objManager,
136 const word& adjointTurbulenceModelName
137)
138:
139 adjointTurbulenceModel
140 (
141 primalVars,
142 adjointVars,
143 objManager,
144 adjointTurbulenceModelName
145 ),
146 IOdictionary
147 (
148 IOobject
149 (
150 "adjointRASProperties",
151 primalVars.U().time().constant(),
152 primalVars.U().db(),
153 IOobject::MUST_READ_IF_MODIFIED,
154 IOobject::NO_WRITE
155 )
156 ),
157
158 objectiveManager_(objManager),
159
160 adjointTurbulence_(get<word>("adjointTurbulence")),
161 printCoeffs_(getOrDefault<Switch>("printCoeffs", false)),
162 coeffDict_(subOrEmptyDict(type + "Coeffs")),
163
164 y_(mesh_),
165
166 adjointTMVariable1Ptr_(nullptr),
167 adjointTMVariable2Ptr_(nullptr),
168 adjointTMVariablesBaseNames_(0),
169 adjointTMVariable1MeanPtr_(nullptr),
170 adjointTMVariable2MeanPtr_(nullptr),
171 adjMomentumBCSourcePtr_( createZeroBoundaryPtr<vector>(mesh_) ),
172 wallShapeSensitivitiesPtr_( createZeroBoundaryPtr<vector>(mesh_) ),
173 wallFloCoSensitivitiesPtr_( createZeroBoundaryPtr<vector>(mesh_) ),
176{}
177
178
179// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
180
182(
183 incompressibleVars& primalVars,
185 objectiveManager& objManager,
186 const word& adjointTurbulenceModelName
187)
188{
189 const IOdictionary dict
190 (
192 (
193 "adjointRASProperties",
194 primalVars.U().time().constant(),
195 primalVars.U().db(),
199 )
200 );
201
202 const word modelType(dict.get<word>("adjointRASModel"));
203
204 Info<< "Selecting adjointRAS turbulence model " << modelType << endl;
205
206 auto* ctorPtr = dictionaryConstructorTable(modelType);
207
208 if (!ctorPtr)
209 {
211 (
212 dict,
213 "adjointRASModel",
214 modelType,
215 *dictionaryConstructorTablePtr_
216 ) << exit(FatalIOError);
217 }
218
220 (
221 ctorPtr
222 (
223 primalVars,
224 adjointVars,
225 objManager,
226 adjointTurbulenceModelName
228 );
229}
230
231
232// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
233
235{
237
238 if (adjointTurbulence_ && mesh_.changing())
239 {
240 y_.correct();
241 }
242}
243
244
246{
247 //if (regIOobject::read())
248
249 // Bit of trickery : we are both IOdictionary ('adjointRASProperties') and
250 // an regIOobject from the adjointTurbulenceModel level. Problem is to
251 // distinguish between the two - we only want to reread the IOdictionary.
252
253 bool ok = IOdictionary::readData
254 (
255 IOdictionary::readStream
256 (
257 IOdictionary::type()
258 )
259 );
261
262 if (ok)
263 {
264 readEntry("adjointTurbulence", adjointTurbulence_);
265
266 if (const dictionary* dictPtr = findDict(type() + "Coeffs"))
267 {
268 coeffDict_ <<= *dictPtr;
269 }
270
271 return true;
272 }
273 else
274 {
275 return false;
276 }
277}
278
279
281{
283 {
284 // if pointer is not set, set it to a zero field
286 (
288 (
290 (
291 "adjointTMVariable1" + type(),
292 mesh_.time().timeName(),
293 mesh_,
296 ),
297 mesh_,
299 )
300 );
301 }
302
303 return adjointTMVariable1Ptr_();
304}
305
306
308{
310 {
311 // if pointer is not set, set it to a zero field
313 (
315 (
317 (
318 "adjointTMVariable2" + type(),
319 mesh_.time().timeName(),
320 mesh_,
323 ),
324 mesh_,
326 )
339 }
340 else
353 }
354 else
355 {
357 }
358}
359
364}
365
370}
371
374{
376}
377
378
380{
381 dimensionSet dims
382 (
384 ? dimViscosity/adjointTMVariable1Ptr_().dimensions()
385 : dimless
386 );
387
388 return
390 (
392 (
393 "nutJacobianTMVar1"+type(),
394 mesh_.time().timeName(),
395 mesh_,
398 ),
399 mesh_,
401 );
402}
403
404
406{
407 dimensionSet dims
408 (
410 ? dimViscosity/adjointTMVariable2Ptr_().dimensions()
411 : dimless
412 );
413
414 return
416 (
418 (
419 "nutJacobianTMVar2"+type(),
420 mesh_.time().timeName(),
421 mesh_,
432(
433 tmp<volScalarField>& dNutdUMult
434) const
435{
436 // Deliberately returning a null pointer in the base class
437 return nullptr;
438}
439
442{
443 return tmp<scalarField>::New(mesh_.boundary()[patchI].size(), Zero);
444}
445
448{
449 return tmp<scalarField>::New(mesh_.boundary()[patchI].size(), Zero);
450}
451
454{
456}
457
458
460{
461 const solverControl& solControl = adjointVars_.getSolverControl();
462 if (solControl.average())
463 {
465 {
467 }
471 }
472 }
473}
474
475
477{
478 const solverControl& solControl = adjointVars_.getSolverControl();
479 if (solControl.doAverageIter())
480 {
481 const label iAverageIter = solControl.averageIter();
482 scalar avIter(iAverageIter);
483 scalar oneOverItP1 = 1./(avIter+1);
484 scalar mult = avIter*oneOverItP1;
486 {
489 + getAdjointTMVariable1Inst()*oneOverItP1;
490 }
492 {
495 + getAdjointTMVariable2Inst()*oneOverItP1;
496 }
497 }
498}
499
500
502{
503 return includeDistance_;
504}
505
506
507// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
508
509} // End namespace incompressible
510} // End namespace Foam
511
512// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
IOdictionary(const IOobject &io, const dictionary *fallback=nullptr)
Construct given an IOobject and optional fallback dictionary content.
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition IOobject.C:450
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition Switch.H:81
const word & constant() const noexcept
Return constant name.
Definition TimePathsI.H:131
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
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
virtual bool readData(Istream &)
The readData function required by regIOobject read operation.
const word & name() const
Name function is needed to disambiguate those inherited from regIOobject and dictionary.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Find and return a sub-dictionary as a copy, otherwise return an empty dictionary.
Definition dictionary.C:521
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.
const dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary pointer if present (and it is a dictionary) otherwise return nullptr...
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,...
dictionary()
Default construct, a top-level empty dictionary.
Definition dictionary.C:68
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
Manages the adjoint mean flow fields and their mean values.
const solverControl & getSolverControl() const
Return const reference to solverControl.
Abstract base class for incompressible turbulence models.
virtual tmp< scalarField > diffusionCoeffVar1(label patchI) const
Diffusion coefficient of the first primal and adjoint turbulence model equation. Needed for some adjo...
virtual tmp< volScalarField > nutJacobianTMVar2() const
Jacobian of nut wrt the second turbulence model variable.
volScalarField & getAdjointTMVariable1Inst()
Return non-constant reference to adjoint turbulence model variable 1.
dictionary coeffDict_
Model coefficients dictionary.
bool includeDistance_
Does the turbulence model include distances and should the adjoint to the distance field be computed.
autoPtr< boundaryVectorField > wallShapeSensitivitiesPtr_
Wall sensitivity term for shape optimisation.
autoPtr< volScalarField > adjointTMVariable1MeanPtr_
Adjoint turbulence model variable 1, mean value.
autoPtr< volScalarField > adjointTMVariable2MeanPtr_
Adjoint turbulence model variable 2, mean value.
volScalarField & getAdjointTMVariable2Inst()
Return non-constant reference to adjoint turbulence model variable 2.
virtual void correct()
Solve the adjoint turbulence equations.
nearWallDist y_
Near wall distance boundary field.
Switch adjointTurbulence_
Turbulence on/off flag.
void restoreInitValues()
Restore field values to the initial ones.
bool changedPrimalSolution_
Has the primal solution changed?
autoPtr< boundaryVectorField > adjMomentumBCSourcePtr_
Source to the adjoint momentum BC emerging from differentiating the turbulence model.
virtual tmp< volScalarField > nutJacobianTMVar1() const
Jacobian of nut wrt the first turbulence model variable.
virtual void printCoeffs()
Print model coefficients.
bool includeDistance() const
Should the adjoint to the eikonal equation be computed.
autoPtr< volScalarField > & getAdjointTMVariable2InstPtr()
Return non-constant autoPtr to adjoint turbulence model variable 2.
autoPtr< volScalarField > adjointTMVariable1Ptr_
Adjoint turbulence model variable 1.
autoPtr< volScalarField > adjointTMVariable2Ptr_
Adjoint turbulence model variable 2.
void computeMeanFields()
Average adjoint fields on the fly.
void resetMeanFields()
Reset mean fields to zero.
void setChangedPrimalSolution()
Set flag of changed primal solution to true.
objectiveManager & objectiveManager_
Reference to the objectiveManager.
volScalarField & getAdjointTMVariable2()
Return non-constant reference to adjoint turbulence model variable 2.
volScalarField & getAdjointTMVariable1()
Return non-constant reference to adjoint turbulence model variable 1.
const wordList & getAdjointTMVariablesBaseNames() const
Return reference to the adjoint turbulence model variables base names.
wordList adjointTMVariablesBaseNames_
Base names of the adjoint fields.
virtual tmp< volVectorField > nutJacobianU(tmp< volScalarField > &dNutdUMult) const
Jacobian of nut wrt the flow velocity.
virtual tmp< scalarField > diffusionCoeffVar2(label patchI) const
Diffusion coefficient of the second primal and adjoint turbulence model equation. Needed for some adj...
autoPtr< boundaryVectorField > wallFloCoSensitivitiesPtr_
Wall sensitivity term for flow control optimisation.
autoPtr< volScalarField > & getAdjointTMVariable1InstPtr()
Return non-constant autoPtr to adjoint turbulence model variable 1.
Switch printCoeffs_
Flag to print the model coeffs at run-time.
virtual bool read()
Read adjointRASProperties dictionary.
static autoPtr< adjointRASModel > New(incompressibleVars &primalVars, incompressibleAdjointMeanFlowVars &adjointVars, objectiveManager &objManager, const word &adjointTurbulenceModelName=adjointTurbulenceModel::typeName)
Return a reference to the selected adjointRAS model.
Abstract base class for incompressible adjoint turbulence models (RAS, LES and laminar).
virtual void correct()=0
Solve the adjoint turbulence equations.
Base class for solution control classes.
const volVectorField & U() const
Return const reference to velocity.
Class for managing objective functions.
constant condensation/saturation model.
void close()
Close Istream.
Base class for solver control classes.
bool doAverageIter() const
Whether or not to add fields of the current iteration to the average fields.
bool storeInitValues() const
Re-initialize.
label & averageIter()
Return average iteration index reference.
bool useAveragedFields() const
Use averaged fields? For solving the adjoint equations or computing sensitivities based on averaged f...
bool average() const
Whether averaging is enabled or not.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
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
U
Definition pEqn.H:72
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition error.H:637
Different types of constants.
Namespace for incompressible adjoint turbulence models.
Namespace for OpenFOAM.
const dimensionSet dimViscosity
List< word > wordList
List of word.
Definition fileName.H:60
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
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< typename GeometricField< Type, fvPatchField, volMesh >::Boundary > createZeroBoundaryPtr(const fvMesh &mesh, bool printAllocation=false)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
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