incompressibleAdjointMeanFlowVars.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-2021 PCOpt/NTUA
9 Copyright (C) 2013-2021 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
31
32// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33
34namespace Foam
35{
36
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
40
41// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
42
44{
48 (
50 mesh_,
51 UaInst(),
52 "phia",
55 );
56
58}
59
61{
62 // Allocate mean fields
63 // Only mean flow here since turbulent quantities
64 // are allocated automatically in RASModelVariables
66 {
67 Info<< "Allocating Mean Adjoint Fields" << endl;
68 paMeanPtr_.reset
69 (
71 (
73 (
74 paInst().name() + "Mean",
76 mesh_,
79 ),
80 paInst()
81 )
82 );
83 UaMeanPtr_.reset
84 (
86 (
88 (
89 UaInst().name() + "Mean",
91 mesh_,
94 ),
95 UaInst()
96 )
97 );
98 phiaMeanPtr_.reset
99 (
101 (
103 (
104 phiaInst().name() + "Mean",
105 mesh_.time().timeName(),
106 mesh_,
109 ),
110 phiaInst()
111 )
112 );
113 }
114}
115
116// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
117
119(
120 fvMesh& mesh,
121 solverControl& SolverControl,
122 incompressibleVars& primalVars
123)
124:
125 variablesSet(mesh, SolverControl.solverDict()),
126 solverControl_(SolverControl),
127 primalVars_(primalVars),
128 paPtr_(nullptr),
129 UaPtr_(nullptr),
130 phiaPtr_(nullptr),
131 paMeanPtr_(nullptr),
132 UaMeanPtr_(nullptr),
133 phiaMeanPtr_(nullptr)
134{
135 setFields();
137}
138
139
140// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
141
143{
144 return primalVars_;
145}
146
148{
150 {
151 return paMeanPtr_();
152 }
153 else
154 {
155 return paPtr_();
156 }
157}
158
159
161{
163 {
164 return paMeanPtr_();
165 }
166 else
167 {
168 return paPtr_();
169 }
170}
171
172
174{
176 {
177 return UaMeanPtr_();
178 }
179 else
180 {
181 return UaPtr_();
182 }
183}
184
185
187{
189 {
190 return UaMeanPtr_();
191 }
192 else
193 {
194 return UaPtr_();
195 }
196}
197
198
200{
202 {
203 return phiaMeanPtr_();
204 }
205 else
206 {
207 return phiaPtr_();
208 }
209}
210
211
213{
215 {
216 return phiaMeanPtr_();
217 }
218 else
219 {
220 return phiaPtr_();
221 }
222}
223
224
226{
227 return solverControl_.average();
228}
229
230
232{
233 return solverControl_;
234}
235
236
238{
242}
243
244
245// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
246
247} // End namespace Foam
248
249// ************************************************************************* //
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Manages the adjoint mean flow fields and their mean values.
solverControl & solverControl_
Reference to the solverControl of the solver allocating the fields.
const volScalarField & pa() const
Return const reference to pressure.
const volVectorField & Ua() const
Return const reference to velocity.
const volVectorField & UaInst() const
Return const reference to velocity.
const surfaceScalarField & phia() const
Return const reference to volume flux.
const solverControl & getSolverControl() const
Return const reference to solverControl.
virtual void nullify()
Nullify all adjoint fields.
void setFields()
Read fields and set turbulence.
const surfaceScalarField & phiaInst() const
Return const reference to volume flux.
bool computeMeanFields() const
Return computeMeanFields bool.
const volScalarField & paInst() const
Return const reference to pressure.
void setMeanFields()
Read mean fields, if necessary.
incompressibleVars & primalVars_
Reference to primal variables.
Base class for solution control classes.
void setFluxRequired(const word &name) const
Get flux-required for given name, or default.
Base class for solver control classes.
Definition: solverControl.H:52
bool useAveragedFields() const
bool average() const
Whether averaging is enabled or not.
Base class for creating a set of variables.
Definition: variablesSet.H:50
bool useSolverNameForFields_
Append the solver name to the variables names?
Definition: variablesSet.H:98
static void setField(autoPtr< GeometricField< Type, fvPatchField, volMesh > > &fieldPtr, const fvMesh &mesh, const word &baseName, const word &solverName, const bool useSolverNameForFields)
Read vol fields.
word solverName_
Solver name owning the variables set.
Definition: variablesSet.H:95
fvMesh & mesh_
Reference to the mesh database.
Definition: variablesSet.H:92
static void setFluxField(autoPtr< surfaceScalarField > &fieldPtr, const fvMesh &mesh, const volVectorField &velocity, const word &baseName, const word &solverName, const bool useSolverNameForFields)
Set flux field.
Definition: variablesSet.C:97
static void nullifyField(GeometricField< Type, PatchField, GeoMesh > &fieldPtr)
Nullify field and old times, if present.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
Namespace for OpenFOAM.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59