incompressibleVars.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-2020 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 "incompressibleVars.H"
31 #include "createZeroField.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
41 
42 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
43 
45 {
49  (
50  phiPtr_,
51  mesh_,
52  UInst(),
53  "phi",
56  );
57 
58  mesh_.setFluxRequired(pPtr_->name());
59 
60  // if required, correct boundary conditions of mean flow fields here in
61  // order to have the correct bcs for e.g. turbulence models that follow.
62  // NOTE: phi correction depends on the solver (includes for instance
63  // Rhie-Chow interpolation). This needs to be implemented within
64  // incompressiblePrimalSolver
66  {
68  }
69 
71  (
73  );
74  turbulence_.reset
75  (
77  (
78  UInst(),
79  phiInst(),
81  ).ptr()
82  );
83  RASModelVariables_.reset
84  (
86  (
87  mesh_,
89  ).ptr()
90  );
93  {
95  }
96 }
97 
98 
100 {
101  // Store init fields
102  // only mean flow here since turbulent quantities
103  // are allocated automatically in RASModelVariables
105  {
106  pInitPtr_.reset(new volScalarField(pInst().name() + "Init", pInst()));
107  UInitPtr_.reset(new volVectorField(UInst().name() + "Init", UInst()));
108  phiInitPtr_.reset
109  (
110  new surfaceScalarField(phiInst().name() + "Init", phiInst())
111  );
112  }
113 }
114 
115 
117 {
118  // Allocate mean fields
119  // only mean flow here since turbulent quantities
120  // are allocated automatically in RASModelVariables
121  if (solverControl_.average())
122  {
123  Info<< "Allocating Mean Primal Fields" << endl;
124  pMeanPtr_.reset
125  (
126  new volScalarField
127  (
128  IOobject
129  (
130  pInst().name()+"Mean",
131  mesh_.time().timeName(),
132  mesh_,
135  ),
136  pInst()
137  )
138  );
139  UMeanPtr_.reset
140  (
141  new volVectorField
142  (
143  IOobject
144  (
145  UInst().name()+"Mean",
146  mesh_.time().timeName(),
147  mesh_,
150  ),
151  UInst()
152  )
153  );
154  phiMeanPtr_.reset
155  (
157  (
158  IOobject
159  (
160  phiInst().name()+"Mean",
161  mesh_.time().timeName(),
162  mesh_,
165  ),
166  phiInst()
167  )
168  );
169 
170  // Correct boundary conditions if necessary
172  {
173  pMeanPtr_().correctBoundaryConditions();
174  UMeanPtr_().correctBoundaryConditions();
175  }
176  }
177 }
178 
179 
181 {
182  // Turbulence model always reads fields with the prescribed name
183  // If a custom name is supplied, check whether this field exists,
184  // copy it to the field known by the turbulence model
185  // and re-name the latter
187  {
188  incompressible::RASModelVariables& rasVars = RASModelVariables_();
189  if (rasVars.hasTMVar1())
190  {
191  renameTurbulenceField(rasVars.TMVar1Inst(), solverName_);
192  }
193  if (rasVars.hasTMVar2())
194  {
195  renameTurbulenceField(rasVars.TMVar2Inst(), solverName_);
196  }
197  if (rasVars.hasNut())
198  {
199  renameTurbulenceField(rasVars.nutRefInst(), solverName_);
200  }
201  }
202 }
203 
204 
206 {
207  Info<< "Correcting (U,p) boundary conditions " << endl;
210  if (solverControl_.average())
211  {
212  pMeanPtr_().correctBoundaryConditions();
213  UMeanPtr_().correctBoundaryConditions();
214  }
215 }
216 
217 
219 {
220  // If required, correct boundary conditions of turbulent fields.
221  // Includes the correction of boundary conditions for averaged fields,
222  // if any
223  Info<< "Correcting boundary conditions of turbulent fields" << endl;
224  RASModelVariables_().correctBoundaryConditions(turbulence_());
225 }
226 
227 
228 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
229 
231 (
232  fvMesh& mesh,
233  solverControl& SolverControl
234 )
235 :
236  variablesSet(mesh, SolverControl.solverDict()),
237  solverControl_(SolverControl),
238  pPtr_(nullptr),
239  UPtr_(nullptr),
240  phiPtr_(nullptr),
241  laminarTransportPtr_(nullptr),
242  turbulence_(nullptr),
243  RASModelVariables_(nullptr),
244 
245  pInitPtr_(nullptr),
246  UInitPtr_(nullptr),
247  phiInitPtr_(nullptr),
248 
249  pMeanPtr_(nullptr),
250  UMeanPtr_(nullptr),
251  phiMeanPtr_(nullptr),
252 
253  correctBoundaryConditions_
254  (
255  SolverControl.solverDict().subOrEmptyDict("fieldReconstruction").
256  getOrDefault<bool>("reconstruct", false)
257  )
258 {
260  setInitFields();
261  setMeanFields();
262 }
263 
264 
266 (
267  const incompressibleVars& vs
268 )
269 :
270  variablesSet(vs.mesh_, vs.solverControl_.solverDict()),
271  solverControl_(vs.solverControl_),
272  pPtr_(allocateRenamedField(vs.pPtr_)),
273  UPtr_(allocateRenamedField(vs.UPtr_)),
274  phiPtr_(allocateRenamedField(vs.phiPtr_)),
275  laminarTransportPtr_(nullptr),
276  turbulence_(nullptr),
277  RASModelVariables_(vs.RASModelVariables_.clone()),
278 
279  pInitPtr_(allocateRenamedField(vs.pInitPtr_)),
280  UInitPtr_(allocateRenamedField(vs.UInitPtr_)),
281  phiInitPtr_(allocateRenamedField(vs.phiInitPtr_)),
282 
283  pMeanPtr_(allocateRenamedField(vs.pMeanPtr_)),
284  UMeanPtr_(allocateRenamedField(UMeanPtr_)),
285  phiMeanPtr_(allocateRenamedField(vs.phiMeanPtr_)),
286 
287  correctBoundaryConditions_(vs.correctBoundaryConditions_)
288 {
289  DebugInfo
290  << "Calling incompressibleVars copy constructor" << endl;
291 }
292 
293 
295 {
296  DebugInfo
297  << "Calling incompressibleVars::clone" << endl;
299  return autoPtr<variablesSet>(new incompressibleVars(*this));
300 }
301 
302 
303 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
304 
306 {
308  {
309  return pMeanPtr_();
310  }
311  else
312  {
313  return pPtr_();
314  }
315 }
316 
317 
319 {
321  {
322  return pMeanPtr_();
323  }
324  else
325  {
326  return pPtr_();
327  }
328 }
329 
330 
332 {
334  {
335  return UMeanPtr_();
336  }
337  else
338  {
339  return UPtr_();
340  }
341 }
342 
343 
345 {
347  {
348  return UMeanPtr_();
349  }
350  else
351  {
352  return UPtr_();
353  }
354 }
355 
356 
358 {
360  {
361  return phiMeanPtr_();
362  }
363  else
364  {
365  return phiPtr_();
366  }
367 }
368 
370 {
372  {
373  return phiMeanPtr_();
374  }
375  else
376  {
377  return phiPtr_();
378  }
379 }
380 
381 
383 {
385  {
386  Info<< "Restoring field values to initial ones" << endl;
387  pInst() == pInitPtr_();
388  UInst() == UInitPtr_();
389  phiInst() == phiInitPtr_();
390  RASModelVariables_().restoreInitValues();
391  }
392 }
393 
394 
396 {
397  if (solverControl_.average())
398  {
399  Info<< "Resetting mean fields to zero" << endl;
400 
401  // Reset fields to zero
402  pMeanPtr_() == dimensionedScalar(pInst().dimensions(), Zero);
403  UMeanPtr_() == dimensionedVector(UInst().dimensions(), Zero);
404  phiMeanPtr_() == dimensionedScalar(phiInst().dimensions(), Zero);
405  RASModelVariables_().resetMeanFields();
407  // Reset averaging iteration index to 0
409  }
410 }
411 
412 
414 {
416  {
417  Info<< "Averaging fields" << endl;
418  label& iAverageIter = solverControl_.averageIter();
419  scalar avIter(iAverageIter);
420  scalar oneOverItP1 = 1./(avIter + 1);
421  scalar mult = avIter*oneOverItP1;
422  pMeanPtr_() == pMeanPtr_()*mult + pInst()*oneOverItP1;
423  UMeanPtr_() == UMeanPtr_()*mult + UInst()*oneOverItP1;
424  phiMeanPtr_() == phiMeanPtr_()*mult + phiInst()*oneOverItP1;
425  RASModelVariables_().computeMeanFields();
426  ++iAverageIter;
427  }
428 }
429 
430 
432 {
434  RASModelVariables_().correctBoundaryConditions(turbulence_());
435 }
436 
439 {
441 }
442 
445 {
446  return solverControl_.average();
447 }
448 
449 
451 {
452  incompressibleVars& incoVars = refCast<incompressibleVars>(vars);
453  // Copy source fields to the ones known by the object
454  swapAndRename(pPtr_, incoVars.pPtr_);
455  swapAndRename(UPtr_, incoVars.UPtr_);
456  swapAndRename(phiPtr_, incoVars.phiPtr_);
458  // Transfer turbulent fields. Copies fields since original fields are
459  // not owned by RASModelVariables but from the turbulence model
460  RASModelVariables_->transfer(incoVars.RASModelVariables()());
461 }
462 
463 
464 bool incompressibleVars::write() const
465 {
466  // Write dummy fields, for continuation only
468  {
469  if (RASModelVariables_().hasTMVar1())
470  {
471  createZeroFieldPtr<scalar>
472  (
473  mesh_,
474  RASModelVariables_().TMVar1BaseName(),
475  RASModelVariables_().TMVar1Inst().dimensions()
476  )().write();
477  }
478  if (RASModelVariables_().hasTMVar2())
479  {
480  createZeroFieldPtr<scalar>
481  (
482  mesh_,
483  RASModelVariables_().TMVar2BaseName(),
484  RASModelVariables_().TMVar2Inst().dimensions()
485  )().write();
486  }
487  if (RASModelVariables_().hasNut())
488  {
489  createZeroFieldPtr<scalar>
490  (
491  mesh_,
492  RASModelVariables_().nutBaseName(),
493  RASModelVariables_().nutRefInst().dimensions()
494  )().write();
495  }
496 
497  return true;
498  }
499 
500  return false;
501 }
502 
503 
504 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
505 
506 } // End namespace Foam
507 
508 // ************************************************************************* //
autoPtr< singlePhaseTransportModel > laminarTransportPtr_
autoPtr< volVectorField > UPtr_
const volScalarField & pInst() const
Return const reference to pressure.
*virtual void transfer(variablesSet &vars)
Transfer the fields of another variablesSet to this.
autoPtr< volScalarField > pMeanPtr_
Manage mean fields. Turbulent mean fields are managed in RASModelVariables.
static autoPtr< IncompressibleTurbulenceModel > New(const volVectorField &U, const surfaceScalarField &phi, const TransportModel &transportModel, const word &propertiesName=turbulenceModel::propertiesName)
Return a reference to the selected turbulence model.
void resetMeanFields()
Reset mean fields to zero.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const surfaceScalarField & phiInst() const
Return const reference to volume flux.
void setMeanFields()
Set mean fields if necessary.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
incompressibleVars(fvMesh &mesh, solverControl &SolverControl)
Construct from mesh.
void correctNonTurbulentBoundaryConditions()
Update boundary conditions of mean-flow.
void setFluxRequired(const word &name) const
Set flux-required for given name (mutable)
void setInitFields()
Set initial fields if necessary.
static autoPtr< RASModelVariables > New(const fvMesh &mesh, const solverControl &SolverControl)
Return a reference to the selected turbulence model.
void correctTurbulentBoundaryConditions()
Update boundary conditions of turbulent fields.
autoPtr< surfaceScalarField > phiInitPtr_
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
autoPtr< volScalarField > pPtr_
Fields involved in the solution of the incompressible NS equations.
virtual autoPtr< variablesSet > clone() const
Clone the incompressibleVars.
const surfaceScalarField & phi() const
Return const reference to volume flux.
Base class for solution control classes.
dynamicFvMesh & mesh
void renameTurbulenceField(GeometricField< Type, fvPatchField, volMesh > &baseField, const word &solverName)
Turbulence model always reads fields with the prescribed name If a custom name is supplied...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const singlePhaseTransportModel & laminarTransport() const
Return const reference to transport model.
autoPtr< surfaceScalarField > phiMeanPtr_
const volScalarField & p() const
Return const reference to pressure.
bool useSolverNameForFields_
Append the solver name to the variables names?
Definition: variablesSet.H:107
Base class for solver control classes.
Definition: solverControl.H:45
void setFields()
Read fields and set turbulence.
Reading is optional [identical to LAZY_READ].
autoPtr< incompressible::RASModelVariables > RASModelVariables_
autoPtr< volVectorField > UMeanPtr_
solverControl & solverControl_
Reference to the solverControl of the solver allocating the fields.
#define DebugInfo
Report an information message using Foam::Info.
const volVectorField & U() const
Return const reference to velocity.
bool storeInitValues() const
Return storeInitValues bool.
void renameTurbulenceFields()
Rename turbulence fields if necessary.
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)
autoPtr< surfaceScalarField > phiPtr_
bool doAverageIter() const
Whether or not to add fields of the current iteration to the average fields.
autoPtr< volVectorField > UInitPtr_
const autoPtr< incompressible::RASModelVariables > & RASModelVariables() const
Return const reference to the turbulence model variables.
void computeMeanFields()
Compute mean fields on the fly.
void restoreInitValues()
Restore field values to the initial ones.
Base class for creating a set of variables.
Definition: variablesSet.H:43
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
bool correctBoundaryConditions_
Update boundary conditions upon construction.
bool average() const
Whether averaging is enabled or not.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Automatically write from objectRegistry::writeObject()
void correctBoundaryConditions()
correct boundaryconditions for all volFields
autoPtr< volScalarField > pInitPtr_
Keep a copy of the initial field values if necessary.
void correctBoundaryConditions()
Correct boundary field.
messageStream Info
Information stream (stdout output on master, null elsewhere)
autoPtr< incompressible::turbulenceModel > turbulence_
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
label & averageIter()
Return average iteration index reference.
bool write() const
Write dummy turbulent fields to allow for continuation in multi-point, turbulent runs.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const volVectorField & UInst() const
Return const reference to velocity.
bool useAveragedFields() const
Use averaged fields? For solving the adjoint equations or computing sensitivities based on averaged f...
A simple single-phase transport model based on viscosityModel.
bool storeInitValues() const
Re-initialize.
void swapAndRename(autoPtr< GeometricField< Type, PatchField, GeoMesh >> &p1, autoPtr< GeometricField< Type, PatchField, GeoMesh >> &p2)
Swap autoPtrs and rename managed fields.
fvMesh & mesh_
Reference to the mesh database.
Definition: variablesSet.H:97
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
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127