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-------------------------------------------------------------------------------
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
37{
38namespace incompressibleAdjoint
39{
40
41// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42
46(
50);
51
52
53// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
54
56{
57 if (printCoeffs_)
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 (
74 (
76 (
79 mesh_,
82 ),
84 )
85 );
86 }
87
89 {
91 (
93 (
95 (
98 mesh_,
101 ),
103 )
104 );
105 }
106 }
107}
108
109
110// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
111
113(
114 const word& type,
115 incompressibleVars& primalVars,
117 objectiveManager& objManager,
118 const word& adjointTurbulenceModelName
119)
120:
122 (
123 primalVars,
124 adjointVars,
125 objManager,
126 adjointTurbulenceModelName
127 ),
129 (
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 (
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 (
270 (
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 (
297 (
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 {
325 }
326}
327
328
329
331{
333 {
335 }
336 else
337 {
339 }
340}
341
342
344{
346}
347
348
350{
352}
353
354
356{
358}
359
360
362{
363 return
365 (
367 (
368 "nutJacobianTMVar1"+type(),
369 mesh_.time().timeName(),
370 mesh_,
373 ),
374 mesh_,
376 (
377 nut().dimensions()/adjointTMVariable1Ptr_().dimensions(),
378 Zero
379 )
380 );
381}
382
383
385{
386 return
388 (
390 (
391 "nutJacobianTMVar2"+type(),
392 mesh_.time().timeName(),
393 mesh_,
396 ),
397 mesh_,
399 (
400 nut().dimensions()/adjointTMVariable2Ptr_().dimensions(),
401 Zero
402 )
403 );
404}
405
406
408{
409 return tmp<scalarField>::New(mesh_.boundary()[patchI].size(), Zero);
410}
411
412
414{
415 return tmp<scalarField>::New(mesh_.boundary()[patchI].size(), Zero);
416}
417
418
420{
422}
423
424
426{
427 const solverControl& solControl = adjointVars_.getSolverControl();
428 if (solControl.average())
429 {
431 {
434 }
436 {
439 }
440 }
441}
442
443
445{
446 const solverControl& solControl = adjointVars_.getSolverControl();
447 if (solControl.doAverageIter())
448 {
449 const label iAverageIter = solControl.averageIter();
450 scalar avIter(iAverageIter);
451 scalar oneOverItP1 = 1./(avIter+1);
452 scalar mult = avIter*oneOverItP1;
454 {
457 + getAdjointTMVariable1Inst()*oneOverItP1;
458 }
460 {
463 + getAdjointTMVariable2Inst()*oneOverItP1;
464 }
465 }
466}
467
468
470{
471 return includeDistance_;
472}
473
474
475// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
476
477} // End namespace incompressible
478} // End namespace Foam
479
480// ************************************************************************* //
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...
Definition: IOdictionary.H:57
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const Time & time() const
Return Time associated with the objectRegistry.
Definition: IOobject.C:506
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:500
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:180
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:78
const word & constant() const
Return constant name.
Definition: TimePathsI.H:96
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
const word & name() const
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX)
Find and return a sub-dictionary pointer if present.
Definition: dictionaryI.H:127
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:712
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
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
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.
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.
bool changedPrimalSolution_
Has the primal solution changed?
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.
volScalarField & getAdjointTMVariable2()
Return non-constant reference to adjoint turbulence model variable 2.
volScalarField & getAdjointTMVariable1()
Return non-constant reference to adjoint turbulence model variable 1.
wordList adjointTMVariablesBaseNames_
Base names of the adjoint fields.
virtual tmp< scalarField > diffusionCoeffVar2(label patchI) const
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.
Abstract base class for incompressible adjoint turbulence models (RAS, LES and laminar).
virtual void correct()=0
Solve the adjoint turbulence equations.
virtual const volScalarField & nut() const
Return the turbulence viscosity.
Base class for solution control classes.
const volVectorField & U() const
Return const reference to velocity.
virtual void correct()
Correct for mesh geom/topo changes.
Definition: nearWallDist.C:113
class for managing incompressible objective functions.
constant condensation/saturation model.
bool changing() const noexcept
Is mesh changing (topology changing and/or moving)
Definition: polyMesh.H:562
virtual bool readData(Istream &)
ReadData function required for regIOobject read operation. Note:
Base class for solver control classes.
Definition: solverControl.H:52
bool doAverageIter() const
label & averageIter()
Return average iteration index reference.
bool useAveragedFields() const
bool average() const
Whether averaging is enabled or not.
A class for managing temporary objects.
Definition: tmp.H:65
type
Volume classification types.
Definition: volumeType.H:66
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimless
Dimensionless.
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: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
IOerror FatalIOError
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
autoPtr< typename GeometricField< Type, fvPatchField, volMesh >::Boundary > createZeroBoundaryPtr(const fvMesh &mesh, bool printAllocation=false)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
static constexpr char close
Definition: FlatOutput.H:74