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,
48  adjointRASModel,
49  adjointTurbulenceModel
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  (
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 :
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(),
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  adjointTMVariable1MeanPtr_(nullptr),
151  adjointTMVariable2MeanPtr_(nullptr),
152  adjMomentumBCSourcePtr_( createZeroBoundaryPtr<vector>(mesh_) ),
153  wallShapeSensitivitiesPtr_( createZeroBoundaryPtr<vector>(mesh_) ),
154  wallFloCoSensitivitiesPtr_( createZeroBoundaryPtr<vector>(mesh_) ),
155  includeDistance_(false),
156  changedPrimalSolution_(true)
157 {}
158 
159 
160 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
161 
163 (
164  incompressibleVars& primalVars,
166  objectiveManager& objManager,
167  const word& adjointTurbulenceModelName
168 )
169 {
170  const IOdictionary dict
171  (
172  IOobject
173  (
174  "adjointRASProperties",
175  primalVars.U().time().constant(),
176  primalVars.U().db(),
179  false // Do not register
180  )
181  );
182 
183  const word modelType(dict.get<word>("adjointRASModel"));
184 
185  Info<< "Selecting adjointRAS turbulence model " << modelType << endl;
186 
187  auto* ctorPtr = dictionaryConstructorTable(modelType);
188 
189  if (!ctorPtr)
190  {
192  (
193  dict,
194  "adjointRASModel",
195  modelType,
196  *dictionaryConstructorTablePtr_
197  ) << exit(FatalIOError);
198  }
199 
201  (
202  ctorPtr
203  (
204  primalVars,
205  adjointVars,
206  objManager,
207  adjointTurbulenceModelName
208  )
209  );
210 }
211 
212 
213 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
214 
216 {
218 
220  {
221  y_.correct();
222  }
223 }
224 
225 
227 {
228  //if (regIOobject::read())
229 
230  // Bit of trickery : we are both IOdictionary ('adjointRASProperties') and
231  // an regIOobject from the adjointTurbulenceModel level. Problem is to
232  // distinguish between the two - we only want to reread the IOdictionary.
233 
234  bool ok = IOdictionary::readData
235  (
236  IOdictionary::readStream
237  (
239  )
240  );
242 
243  if (ok)
244  {
245  readEntry("adjointTurbulence", adjointTurbulence_);
246 
247  if (const dictionary* dictPtr = findDict(type() + "Coeffs"))
248  {
249  coeffDict_ <<= *dictPtr;
250  }
251 
252  return true;
253  }
254  else
255  {
256  return false;
257  }
258 }
259 
260 
262 {
264  {
265  // if pointer is not set, set it to a zero field
267  (
268  new volScalarField
269  (
270  IOobject
271  (
272  "adjointTMVariable1" + type(),
273  mesh_.time().timeName(),
274  mesh_,
277  ),
278  mesh_,
280  )
281  );
282  }
283 
284  return adjointTMVariable1Ptr_();
285 }
286 
287 
289 {
291  {
292  // if pointer is not set, set it to a zero field
294  (
295  new volScalarField
296  (
297  IOobject
298  (
299  "adjointTMVariable2" + type(),
300  mesh_.time().timeName(),
301  mesh_,
304  ),
305  mesh_,
307  )
308  );
309  }
310 
311  return adjointTMVariable2Ptr_();
312 }
313 
314 
316 {
318  {
320  }
321  else
322  {
323  return getAdjointTMVariable1Inst();
324  }
325 }
326 
327 
328 
330 {
332  {
334  }
335  else
336  {
337  return getAdjointTMVariable2Inst();
338  }
339 }
340 
341 
343 {
344  return adjointTMVariable1Ptr_;
345 }
346 
347 
349 {
350  return adjointTMVariable2Ptr_;
351 }
352 
353 
355 {
356  return
358  (
359  IOobject
360  (
361  "nutJacobianTMVar1"+type(),
362  mesh_.time().timeName(),
363  mesh_,
366  ),
367  mesh_,
369  (
370  nut().dimensions()/adjointTMVariable1Ptr_().dimensions(),
371  Zero
372  )
373  );
374 }
375 
376 
378 {
379  return
381  (
382  IOobject
383  (
384  "nutJacobianTMVar2"+type(),
385  mesh_.time().timeName(),
386  mesh_,
389  ),
390  mesh_,
392  (
393  nut().dimensions()/adjointTMVariable2Ptr_().dimensions(),
394  Zero
395  )
396  );
397 }
398 
399 
401 {
402  return tmp<scalarField>::New(mesh_.boundary()[patchI].size(), Zero);
403 }
404 
405 
407 {
408  return tmp<scalarField>::New(mesh_.boundary()[patchI].size(), Zero);
409 }
410 
411 
413 {
414  changedPrimalSolution_ = true;
415 }
416 
417 
419 {
420  const solverControl& solControl = adjointVars_.getSolverControl();
421  if (solControl.average())
422  {
424  {
427  }
429  {
432  }
433  }
434 }
435 
436 
438 {
439  const solverControl& solControl = adjointVars_.getSolverControl();
440  if (solControl.doAverageIter())
441  {
442  const label iAverageIter = solControl.averageIter();
443  scalar avIter(iAverageIter);
444  scalar oneOverItP1 = 1./(avIter+1);
445  scalar mult = avIter*oneOverItP1;
447  {
450  + getAdjointTMVariable1Inst()*oneOverItP1;
451  }
453  {
456  + getAdjointTMVariable2Inst()*oneOverItP1;
457  }
458  }
459 }
460 
461 
463 {
464  return includeDistance_;
465 }
466 
467 
468 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
469 
470 } // End namespace incompressible
471 } // End namespace Foam
472 
473 // ************************************************************************* //
Foam::dictionary::findDict
dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX)
Find and return a sub-dictionary pointer if present.
Definition: dictionaryI.H:127
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::incompressibleAdjointMeanFlowVars
Manages the adjoint mean flow fields and their mean values.
Definition: incompressibleAdjointMeanFlowVars.H:51
Foam::objectiveManager
class for managing incompressible objective functions.
Definition: objectiveManager.H:54
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::baseIOdictionary::readData
virtual bool readData(Istream &)
The readData function required by regIOobject read operation.
Definition: baseIOdictionary.C:91
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
Foam::solverControl::doAverageIter
bool doAverageIter() const
Definition: solverControlI.H:81
Foam::incompressibleAdjoint::adjointRASModel::getAdjointTMVariable1
volScalarField & getAdjointTMVariable1()
Return non-constant reference to adjoint turbulence model variable 1.
Definition: adjointRASModel.C:315
Foam::incompressibleAdjoint::adjointRASModel::nutJacobianTMVar2
virtual tmp< volScalarField > nutJacobianTMVar2() const
Jacobian of nut wrt the second turbulence model variable.
Definition: adjointRASModel.C:377
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
adjointRASModel.H
Foam::incompressibleAdjoint::addToRunTimeSelectionTable
addToRunTimeSelectionTable(adjointTurbulenceModel, adjointRASModel, adjointTurbulenceModel)
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::incompressibleAdjoint::adjointRASModel::nutJacobianTMVar1
virtual tmp< volScalarField > nutJacobianTMVar1() const
Jacobian of nut wrt the first turbulence model variable.
Definition: adjointRASModel.C:354
Foam::incompressibleAdjoint::adjointRASModel::adjointTMVariable2MeanPtr_
autoPtr< volScalarField > adjointTMVariable2MeanPtr_
Adjoint turbulence model variable 2, mean value.
Definition: adjointRASModel.H:122
Foam::incompressibleAdjoint::adjointRASModel::coeffDict_
dictionary coeffDict_
Model coefficients dictionary.
Definition: adjointRASModel.H:107
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
wallFvPatch.H
Foam::incompressibleAdjoint::adjointRASModel::correct
virtual void correct()
Solve the adjoint turbulence equations.
Definition: adjointRASModel.C:215
Foam::FatalIOError
IOerror FatalIOError
Foam::solverControl
Base class for solver control classes.
Definition: solverControl.H:51
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::incompressibleAdjoint::adjointRASModel::changedPrimalSolution_
bool changedPrimalSolution_
Has the primal solution changed?
Definition: adjointRASModel.H:139
Foam::incompressibleAdjoint::adjointRASModel::adjointTMVariable1MeanPtr_
autoPtr< volScalarField > adjointTMVariable1MeanPtr_
Adjoint turbulence model variable 1, mean value.
Definition: adjointRASModel.H:119
Foam::incompressibleAdjoint::adjointRASModel::diffusionCoeffVar2
virtual tmp< scalarField > diffusionCoeffVar2(label patchI) const
Definition: adjointRASModel.C:406
Foam::baseIOdictionary::name
const word & name() const
Definition: baseIOdictionary.C:85
Foam::incompressibleAdjoint::adjointRASModel::y_
nearWallDist y_
Near wall distance boundary field.
Definition: adjointRASModel.H:110
Foam::incompressibleAdjoint::adjointRASModel::read
virtual bool read()
Read adjointRASProperties dictionary.
Definition: adjointRASModel.C:226
Foam::nearWallDist::correct
virtual void correct()
Correct for mesh geom/topo changes.
Definition: nearWallDist.C:113
Foam::incompressibleAdjoint::adjointRASModel::adjointTMVariable2Ptr_
autoPtr< volScalarField > adjointTMVariable2Ptr_
Adjoint turbulence model variable 2.
Definition: adjointRASModel.H:116
Foam::incompressibleAdjoint::adjointTurbulenceModel::mesh_
const fvMesh & mesh_
Definition: adjointTurbulenceModel.H:93
Foam::incompressibleAdjoint::adjointRASModel::getAdjointTMVariable2Inst
volScalarField & getAdjointTMVariable2Inst()
Return non-constant reference to adjoint turbulence model variable 2.
Definition: adjointRASModel.C:288
Foam::incompressibleAdjoint::defineTypeNameAndDebug
defineTypeNameAndDebug(adjointRASModel, 0)
FatalIOErrorInLookup
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
Foam::incompressibleAdjoint::adjointRASModel::getAdjointTMVariable2
volScalarField & getAdjointTMVariable2()
Return non-constant reference to adjoint turbulence model variable 2.
Definition: adjointRASModel.C:329
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:187
Foam::solverControl::useAveragedFields
bool useAveragedFields() const
Definition: solverControlI.H:94
Foam::incompressibleAdjoint::adjointRASModel::getAdjointTMVariable2InstPtr
autoPtr< volScalarField > & getAdjointTMVariable2InstPtr()
Return non-constant autoPtr to adjoint turbulence model variable 2.
Definition: adjointRASModel.C:348
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:302
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::incompressibleVars::U
const volVectorField & U() const
Return const reference to velocity.
Definition: incompressibleVars.C:331
Foam::incompressibleAdjoint::adjointTurbulenceModel
Abstract base class for incompressible adjoint turbulence models (RAS, LES and laminar).
Definition: adjointTurbulenceModel.H:71
Foam::incompressibleAdjoint::adjointRASModel::setChangedPrimalSolution
void setChangedPrimalSolution()
Set flag of changed primal solution to true.
Definition: adjointRASModel.C:412
Foam::solverControl::averageIter
label & averageIter()
Return average iteration index reference.
Definition: solverControlI.H:63
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::incompressibleAdjoint::adjointRASModel::adjointTurbulence_
Switch adjointTurbulence_
Turbulence on/off flag.
Definition: adjointRASModel.H:101
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::incompressibleAdjoint::adjointRASModel::includeDistance
bool includeDistance() const
Should the adjoint to the eikonal equation be computed.
Definition: adjointRASModel.C:462
Foam::incompressibleAdjoint::adjointTurbulenceModel::correct
virtual void correct()=0
Solve the adjoint turbulence equations.
Definition: adjointTurbulenceModel.C:129
Foam::dimensioned< scalar >
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::incompressibleAdjointMeanFlowVars::getSolverControl
const solverControl & getSolverControl() const
Return const reference to solverControl.
Definition: incompressibleAdjointMeanFlowVars.C:267
Foam::incompressibleAdjoint::defineRunTimeSelectionTable
defineRunTimeSelectionTable(adjointRASModel, dictionary)
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::incompressibleAdjoint::adjointRASModel::getAdjointTMVariable1Inst
volScalarField & getAdjointTMVariable1Inst()
Return non-constant reference to adjoint turbulence model variable 1.
Definition: adjointRASModel.C:261
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:685
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::polyMesh::changing
bool changing() const noexcept
Is mesh changing (topology changing and/or moving)
Definition: polyMesh.H:550
Foam::regIOobject::close
void close()
Close Istream.
Definition: regIOobjectRead.C:171
Foam::incompressibleAdjoint::adjointRASModel::adjointTMVariable1Ptr_
autoPtr< volScalarField > adjointTMVariable1Ptr_
Adjoint turbulence model variable 1.
Definition: adjointRASModel.H:113
Foam::incompressibleAdjoint::adjointRASModel::printCoeffs
virtual void printCoeffs()
Print model coefficients.
Definition: adjointRASModel.C:55
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::incompressibleAdjoint::adjointRASModel::getAdjointTMVariable1InstPtr
autoPtr< volScalarField > & getAdjointTMVariable1InstPtr()
Return non-constant autoPtr to adjoint turbulence model variable 1.
Definition: adjointRASModel.C:342
Foam::incompressibleAdjoint::adjointRASModel::diffusionCoeffVar1
virtual tmp< scalarField > diffusionCoeffVar1(label patchI) const
Definition: adjointRASModel.C:400
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::incompressibleAdjoint::adjointTurbulenceModel::adjointVars_
incompressibleAdjointMeanFlowVars & adjointVars_
Definition: adjointTurbulenceModel.H:91
Foam::IOobject::MUST_READ_IF_MODIFIED
Definition: IOobject.H:186
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
Foam::incompressibleAdjoint::adjointRASModel::resetMeanFields
void resetMeanFields()
Reset mean fields to zero.
Definition: adjointRASModel.C:418
Foam::incompressibleAdjoint::adjointRASModel::includeDistance_
bool includeDistance_
Definition: adjointRASModel.H:136
Foam::incompressibleAdjoint::adjointRASModel::setMeanFields
void setMeanFields()
Set mean fields.
Definition: adjointRASModel.C:64
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::incompressibleAdjoint::adjointTurbulenceModel::nut
virtual const volScalarField & nut() const
Return the turbulence viscosity.
Definition: adjointTurbulenceModel.H:161
Foam::incompressibleAdjoint::adjointRASModel::printCoeffs_
Switch printCoeffs_
Flag to print the model coeffs at run-time.
Definition: adjointRASModel.H:104
Foam::solverControl::average
bool average() const
Whether averaging is enabled or not.
Definition: solverControlI.H:107
Foam::incompressibleVars
Base class for solution control classes.
Definition: incompressibleVars.H:54
Foam::incompressibleAdjoint::adjointRASModel::computeMeanFields
void computeMeanFields()
Average adjoint fields on the fly.
Definition: adjointRASModel.C:437
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
Foam::incompressibleAdjoint::adjointRASModel::New
static autoPtr< adjointRASModel > New(incompressibleVars &primalVars, incompressibleAdjointMeanFlowVars &adjointVars, objectiveManager &objManager, const word &adjointTurbulenceModelName=adjointTurbulenceModel::typeName)
Return a reference to the selected adjointRAS model.
Definition: adjointRASModel.C:163