adjointLaminar.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-2023 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 "adjointLaminar.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace incompressibleAdjoint
38 {
39 namespace adjointRASModels
40 {
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
49 adjointLaminar::adjointLaminar
50 (
51  incompressibleVars& primalVars,
53  objectiveManager& objManager,
54  const word& adjointTurbulenceModelName,
55  const word& modelName
56 )
57 :
59  (
60  modelName,
61  primalVars,
62  adjointVars,
63  objManager,
64  adjointTurbulenceModelName
65  )
66 {}
67 
68 
69 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
70 
72 {
73  const volVectorField& Ua = adjointVars_.Ua();
74  return devReff(Ua);
75 }
76 
77 
79 (
80  const volVectorField& U
81 ) const
82 {
84  (
86  (
87  IOobject
88  (
89  "devRhoReff",
91  mesh_,
94  ),
96  )
97  );
98 }
99 
100 
102 {
103  return
104  (
105  - fvm::laplacian(nuEff(), U)
106  - fvc::div(nuEff()*dev(T(fvc::grad(U))))
107  );
108 }
109 
110 
112 {
114  (
115  IOobject
116  (
117  "adjointMeanFlowSource",
118  runTime_.timeName(),
119  mesh_,
122  ),
123  mesh_,
125  (
126  dimensionSet(0, 1, -2, 0, 0),
127  Zero
128  )
129  );
130 }
131 
132 
134 {
135  // zero contribution
136  return adjMomentumBCSourcePtr_();
137 }
138 
141 {
143 }
144 
147 {
149 }
150 
151 
153 {
155  (
156  IOobject
157  (
158  "adjointEikonalSource" + type(),
159  runTime_.timeName(),
160  mesh_,
163  ),
164  mesh_,
166  );
167 }
168 
169 
171 {
173  (
174  IOobject
175  (
176  "volumeSensTerm" + type(),
177  runTime_.timeName(),
178  mesh_,
181  ),
183  dimensionedTensor(dimensionSet(0, 2, -3, 0, 0), Zero)
184  );
185 }
186 
187 
189 (
190  const word& designVarsName
191 ) const
192 {
194 }
195 
198 {
199  // Does nothing. No fields to nullify
200 }
201 
204 {
205  return adjointRASModel::read();
206 }
207 
208 
210 {
212 }
213 
214 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215 
216 } // End namespace adjointRASModels
217 } // End namespace incompressibleAdjoint
218 } // End namespace Foam
219 
220 // ************************************************************************* //
virtual void correct()
Correct the primal viscosity field. Redundant?
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:85
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
virtual tmp< volScalarField > distanceSensitivities()
Returns zero field.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
Class for managing objective functions.
virtual bool read()
Read adjointRASProperties dictionary.
virtual tmp< fvVectorMatrix > divDevReff(volVectorField &U) const
Return the diffusion term for the momentum equation.
autoPtr< boundaryVectorField > adjMomentumBCSourcePtr_
Source to the adjoint momentum BC emerging from differentiating the turbulence model.
virtual void correct()=0
Solve the adjoint turbulence equations.
Ignore writing from objectRegistry::writeObject()
autoPtr< boundaryVectorField > wallShapeSensitivitiesPtr_
Wall sensitivity term for shape optimisation.
Macros for easy insertion into run-time selection tables.
tmp< volScalarField > nu() const
Return the laminar viscosity.
virtual bool read()
Read adjointRASProperties dictionary.
autoPtr< boundaryVectorField > wallFloCoSensitivitiesPtr_
Wall sensitivity term for flow control optimisation.
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
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Dummy turbulence model for a laminar incompressible flow. Can also be used when the "frozen turbulenc...
virtual void nullify()
Nullify all adjoint turbulence model fields and their old times.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Abstract base class for incompressible turbulence models.
const volVectorField & Ua() const
Return const reference to velocity.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
virtual tmp< volTensorField > FISensitivityTerm()
Returns zero field.
addToRunTimeSelectionTable(adjointRASModel, adjointkOmegaSST, dictionary)
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
virtual const boundaryVectorField & wallFloCoSensitivities()
Returns zero field.
virtual tmp< volSymmTensorField > devReff() const
Return the effective stress tensor, i.e. the adjointLaminar stress.
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
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
virtual const boundaryVectorField & wallShapeSensitivities()
Returns zero field.
Manages the adjoint mean flow fields and their mean values.
dimensionedScalar pow3(const dimensionedScalar &ds)
U
Definition: pEqn.H:72
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label nCells() const noexcept
Number of mesh cells.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
Nothing to be read.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
virtual const boundaryVectorField & adjointMomentumBCSource() const
Returns zero field.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Namespace for OpenFOAM.
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
Definition: SymmTensorI.H:491
virtual tmp< volVectorField > adjointMeanFlowSource()
Source terms to the adjoint momentum equation due to the differentiation of the turbulence model...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
virtual tmp< scalarField > topologySensitivities(const word &designVarsName) const
Returns zero field.