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 -------------------------------------------------------------------------------
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 
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
40 
41 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
42 
44 {
48  (
49  phiaPtr_,
50  mesh_,
51  UaInst(),
52  "phia",
55  );
56 
57  mesh_.setFluxRequired(paPtr_->name());
58 }
59 
60 
62 {
63  // Allocate mean fields
64  // Only mean flow here since turbulent quantities
65  // are allocated automatically in RASModelVariables
66  if (solverControl_.average())
67  {
68  Info<< "Allocating Mean Adjoint Fields" << endl;
69  paMeanPtr_.reset
70  (
71  new volScalarField
72  (
73  IOobject
74  (
75  paInst().name() + "Mean",
76  mesh_.time().timeName(),
77  mesh_,
80  ),
81  paInst()
82  )
83  );
84  UaMeanPtr_.reset
85  (
86  new volVectorField
87  (
88  IOobject
89  (
90  UaInst().name() + "Mean",
91  mesh_.time().timeName(),
92  mesh_,
95  ),
96  UaInst()
97  )
98  );
99  phiaMeanPtr_.reset
100  (
102  (
103  IOobject
104  (
105  phiaInst().name() + "Mean",
106  mesh_.time().timeName(),
107  mesh_,
110  ),
111  phiaInst()
112  )
113  );
114  }
115 }
116 
117 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
118 
119 incompressibleAdjointMeanFlowVars::incompressibleAdjointMeanFlowVars
120 (
121  fvMesh& mesh,
122  solverControl& SolverControl,
123  incompressibleVars& primalVars
124 )
125 :
126  variablesSet(mesh, SolverControl.solverDict()),
127  solverControl_(SolverControl),
128  primalVars_(primalVars),
129  paPtr_(nullptr),
130  UaPtr_(nullptr),
131  phiaPtr_(nullptr),
132  paMeanPtr_(nullptr),
133  UaMeanPtr_(nullptr),
134  phiaMeanPtr_(nullptr)
135 {
137  setMeanFields();
138 }
139 
140 
141 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
142 
144 {
145  return primalVars_;
146 }
147 
149 {
151  {
152  return paMeanPtr_();
153  }
154  else
155  {
156  return paPtr_();
157  }
158 }
159 
160 
162 {
164  {
165  return paMeanPtr_();
166  }
167  else
168  {
169  return paPtr_();
170  }
171 }
172 
173 
175 {
177  {
178  return UaMeanPtr_();
179  }
180  else
181  {
182  return UaPtr_();
183  }
184 }
185 
186 
188 {
190  {
191  return UaMeanPtr_();
192  }
193  else
194  {
195  return UaPtr_();
196  }
197 }
198 
199 
201 {
203  {
204  return phiaMeanPtr_();
205  }
206  else
207  {
208  return phiaPtr_();
209  }
210 }
211 
212 
214 {
216  {
217  return phiaMeanPtr_();
218  }
219  else
220  {
221  return phiaPtr_();
222  }
223 }
224 
227 {
228  return solverControl_.average();
229 }
230 
233 {
234  return solverControl_;
235 }
236 
237 
239 {
243 }
244 
245 
246 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
247 
248 } // End namespace Foam
249 
250 // ************************************************************************* //
void setMeanFields()
Read mean fields, if necessary.
virtual void nullify()
Nullify all adjoint fields.
const volScalarField & paInst() const
Return const reference to pressure.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
autoPtr< volScalarField > paPtr_
Fields involved in the solution of the incompressible adjoint NS equations.
static void nullifyField(GeometricField< Type, PatchField, GeoMesh > &fieldPtr)
Nullify field and old times, if present.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
void setFluxRequired(const word &name) const
Set flux-required for given name (mutable)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
Base class for solution control classes.
incompressibleVars & primalVars_
Reference to primal variables.
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const volScalarField & pa() const
Return const reference to pressure.
const volVectorField & Ua() const
Return const reference to velocity.
bool useSolverNameForFields_
Append the solver name to the variables names?
Definition: variablesSet.H:107
Base class for solver control classes.
Definition: solverControl.H:45
Reading is optional [identical to LAZY_READ].
autoPtr< volScalarField > paMeanPtr_
Mean Adjoint Fields. Actual averaging is done in the incompressibleAdjointVars class to take care of ...
bool computeMeanFields() const
Return computeMeanFields bool.
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
defineTypeNameAndDebug(combustionModel, 0)
const solverControl & getSolverControl() const
Return const reference to solverControl.
Manages the adjoint mean flow fields and their mean values.
const volVectorField & UaInst() const
Return const reference to velocity.
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:90
solverControl & solverControl_
Reference to the solverControl of the solver allocating the fields.
const surfaceScalarField & phiaInst() const
Return const reference to volume flux.
bool average() const
Whether averaging is enabled or not.
Automatically write from objectRegistry::writeObject()
messageStream Info
Information stream (stdout output on master, null elsewhere)
const surfaceScalarField & phia() const
Return const reference to volume flux.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
bool useAveragedFields() const
Use averaged fields? For solving the adjoint equations or computing sensitivities based on averaged f...
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
fvMesh & mesh_
Reference to the mesh database.
Definition: variablesSet.H:97
void setFields()
Read fields and set turbulence.
static void setField(autoPtr< GeometricField< Type, fvPatchField, volMesh >> &fieldPtr, const fvMesh &mesh, const word &baseName, const word &solverName, const bool useSolverNameForFields)
Read vol fields.
Namespace for OpenFOAM.
word solverName_
Solver name owning the variables set.
Definition: variablesSet.H:102