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  (
85  "devRhoReff",
87  (
89  )
90  );
91 }
92 
93 
95 {
96  return
97  (
98  - fvm::laplacian(nuEff(), U)
99  - fvc::div(nuEff()*dev(T(fvc::grad(U))))
100  );
101 }
102 
103 
105 {
106  return volVectorField::New
107  (
108  "adjointMeanFlowSource",
110  mesh_,
111  dimensionedVector(dimensionSet(0, 1, -2, 0, 0), Zero)
112  );
113 }
114 
115 
117 {
118  // zero contribution
119  return adjMomentumBCSourcePtr_();
120 }
121 
124 {
126 }
127 
130 {
132 }
133 
134 
136 {
137  return volScalarField::New
138  (
139  "adjointEikonalSource" + type(),
141  mesh_,
143  );
144 }
145 
146 
148 {
149  return volTensorField::New
150  (
151  "volumeSensTerm" + type(),
154  dimensionedTensor(dimensionSet(0, 2, -3, 0, 0), Zero)
155  );
156 }
157 
158 
160 (
161  const word& designVarsName
162 ) const
163 {
165 }
166 
169 {
170  // Does nothing. No fields to nullify
171 }
172 
175 {
176  return adjointRASModel::read();
177 }
178 
179 
181 {
183 }
184 
185 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
186 
187 } // End namespace adjointRASModels
188 } // End namespace incompressibleAdjoint
189 } // End namespace Foam
190 
191 // ************************************************************************* //
virtual void correct()
Correct the primal viscosity field. Redundant?
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.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
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.
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.
static tmp< GeometricField< symmTensor, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< symmTensor >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
virtual tmp< volSymmTensorField > devReff() const
Return the effective stress tensor, i.e. the adjointLaminar stress.
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
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.
Do not request registration (bool: false)
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.