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-2019 PCOpt/NTUA
9 Copyright (C) 2013-2019 FOSS GP
10 Copyright (C) 2019 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 "adjointLaminar.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace incompressibleAdjoint
38{
39namespace adjointRASModels
40{
41
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43
46
47// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48
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 (
88 (
89 "devRhoReff",
91 mesh_,
94 ),
96 )
97 );
98}
99
100
102{
103 return
104 (
106 - fvc::div(nuEff()*dev(T(fvc::grad(U))))
107 );
108}
109
110
112{
114 (
116 (
117 "adjointMeanFlowSource",
119 mesh_,
122 ),
123 mesh_,
125 (
126 dimensionSet(0, 1, -2, 0, 0),
127 Zero
128 )
129 );
130}
131
132
134{
135 // zero contribution
137}
138
139
141{
143}
144
145
147{
149}
150
151
153{
155 (
157 (
158 "adjointEikonalSource" + type(),
160 mesh_,
163 ),
164 mesh_,
166 );
167}
168
169
171{
173 (
175 (
176 "volumeSensTerm" + type(),
178 mesh_,
181 ),
182 mesh_,
183 dimensionedTensor(dimensionSet(0, 2, -3, 0, 0), Zero)
184 );
185}
186
187
189{
190 // Does nothing. No fields to nullify
191}
192
193
195{
196 return adjointRASModel::read();
197}
198
199
201{
203}
204
205// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
206
207} // End namespace adjointRASModels
208} // End namespace incompressibleAdjoint
209} // End namespace Foam
210
211// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Generic GeometricBoundaryField class.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
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
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
Manages the adjoint mean flow fields and their mean values.
const volVectorField & Ua() const
Return const reference to velocity.
Abstract base class for incompressible turbulence models.
autoPtr< boundaryVectorField > wallShapeSensitivitiesPtr_
Wall sensitivity term for shape optimisation.
autoPtr< boundaryVectorField > adjMomentumBCSourcePtr_
autoPtr< boundaryVectorField > wallFloCoSensitivitiesPtr_
Wall sensitivity term for flow control optimisation.
virtual bool read()
Read adjointRASProperties dictionary.
Dummy turbulence model for a laminar incompressible flow. Can also be used when the "frozen turbulenc...
virtual tmp< fvVectorMatrix > divDevReff(volVectorField &U) const
Return the diffusion term for the momentum equation.
virtual tmp< volTensorField > FISensitivityTerm()
Returns zero field.
virtual void correct()
Correct the primal viscosity field. Redundant?
virtual const boundaryVectorField & adjointMomentumBCSource() const
Returns zero field.
virtual tmp< volSymmTensorField > devReff() const
Return the effective stress tensor, i.e. the adjointLaminar stress.
virtual void nullify()
Nullify all adjoint turbulence model fields and their old times.
virtual const boundaryVectorField & wallShapeSensitivities()
Returns zero field.
virtual const boundaryVectorField & wallFloCoSensitivities()
Returns zero field.
virtual tmp< volScalarField > distanceSensitivities()
Returns zero field.
virtual bool read()
Read adjointRASProperties dictionary.
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
tmp< volScalarField > nu() const
Return the laminar viscosity.
virtual void correct()=0
Solve the adjoint turbulence equations.
Base class for solution control classes.
class for managing incompressible objective functions.
A class for managing temporary objects.
Definition: tmp.H:65
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
const volScalarField & T
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
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
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131