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_,
136  ),
137  pInst()
138  )
139  );
140  UMeanPtr_.reset
141  (
142  new volVectorField
143  (
144  IOobject
145  (
146  UInst().name()+"Mean",
147  mesh_.time().timeName(),
148  mesh_,
152  ),
153  UInst()
154  )
155  );
156  phiMeanPtr_.reset
157  (
159  (
160  IOobject
161  (
162  phiInst().name()+"Mean",
163  mesh_.time().timeName(),
164  mesh_,
168  ),
169  phiInst()
170  )
171  );
172 
173  // Correct boundary conditions if necessary
175  {
176  pMeanPtr_().correctBoundaryConditions();
177  UMeanPtr_().correctBoundaryConditions();
178  }
179  }
180 }
181 
182 
184 {
185  // Turbulence model always reads fields with the prescribed name
186  // If a custom name is supplied, check whether this field exists,
187  // copy it to the field known by the turbulence model
188  // and re-name the latter
190  {
191  incompressible::RASModelVariables& rasVars = RASModelVariables_();
192  if (rasVars.hasTMVar1())
193  {
194  renameTurbulenceField(rasVars.TMVar1Inst(), solverName_);
195  }
196  if (rasVars.hasTMVar2())
197  {
198  renameTurbulenceField(rasVars.TMVar2Inst(), solverName_);
199  }
200  if (rasVars.hasNut())
201  {
202  renameTurbulenceField(rasVars.nutRefInst(), solverName_);
203  }
204  }
205 }
206 
207 
209 {
210  Info<< "Correcting (U,p) boundary conditions " << endl;
213  if (solverControl_.average())
214  {
215  pMeanPtr_().correctBoundaryConditions();
216  UMeanPtr_().correctBoundaryConditions();
217  }
218 }
219 
220 
222 {
223  // If required, correct boundary conditions of turbulent fields.
224  // Includes the correction of boundary conditions for averaged fields,
225  // if any
226  Info<< "Correcting boundary conditions of turbulent fields" << endl;
227  RASModelVariables_().correctBoundaryConditions(turbulence_());
228 }
229 
230 
231 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
232 
234 (
235  fvMesh& mesh,
236  solverControl& SolverControl
237 )
238 :
239  variablesSet(mesh, SolverControl.solverDict()),
240  solverControl_(SolverControl),
241  pPtr_(nullptr),
242  UPtr_(nullptr),
243  phiPtr_(nullptr),
244  laminarTransportPtr_(nullptr),
245  turbulence_(nullptr),
246  RASModelVariables_(nullptr),
247 
248  pInitPtr_(nullptr),
249  UInitPtr_(nullptr),
250  phiInitPtr_(nullptr),
251 
252  pMeanPtr_(nullptr),
253  UMeanPtr_(nullptr),
254  phiMeanPtr_(nullptr),
255 
256  correctBoundaryConditions_
257  (
258  SolverControl.solverDict().subOrEmptyDict("fieldReconstruction").
259  getOrDefault<bool>("reconstruct", false)
260  )
261 {
263  setInitFields();
264  setMeanFields();
265 }
266 
267 
269 (
270  const incompressibleVars& vs
271 )
272 :
273  variablesSet(vs.mesh_, vs.solverControl_.solverDict()),
274  solverControl_(vs.solverControl_),
275  pPtr_(allocateRenamedField(vs.pPtr_)),
276  UPtr_(allocateRenamedField(vs.UPtr_)),
277  phiPtr_(allocateRenamedField(vs.phiPtr_)),
278  laminarTransportPtr_(nullptr),
279  turbulence_(nullptr),
280  RASModelVariables_(vs.RASModelVariables_.clone()),
281 
282  pInitPtr_(allocateRenamedField(vs.pInitPtr_)),
283  UInitPtr_(allocateRenamedField(vs.UInitPtr_)),
284  phiInitPtr_(allocateRenamedField(vs.phiInitPtr_)),
285 
286  pMeanPtr_(allocateRenamedField(vs.pMeanPtr_)),
287  UMeanPtr_(allocateRenamedField(UMeanPtr_)),
288  phiMeanPtr_(allocateRenamedField(vs.phiMeanPtr_)),
289 
290  correctBoundaryConditions_(vs.correctBoundaryConditions_)
291 {
292  DebugInfo
293  << "Calling incompressibleVars copy constructor" << endl;
294 }
295 
296 
298 {
299  DebugInfo
300  << "Calling incompressibleVars::clone" << endl;
302  return autoPtr<variablesSet>(new incompressibleVars(*this));
303 }
304 
305 
306 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
307 
309 {
311  {
312  return pMeanPtr_();
313  }
314  else
315  {
316  return pPtr_();
317  }
318 }
319 
320 
322 {
324  {
325  return pMeanPtr_();
326  }
327  else
328  {
329  return pPtr_();
330  }
331 }
332 
333 
335 {
337  {
338  return UMeanPtr_();
339  }
340  else
341  {
342  return UPtr_();
343  }
344 }
345 
346 
348 {
350  {
351  return UMeanPtr_();
352  }
353  else
354  {
355  return UPtr_();
356  }
357 }
358 
359 
361 {
363  {
364  return phiMeanPtr_();
365  }
366  else
367  {
368  return phiPtr_();
369  }
370 }
371 
373 {
375  {
376  return phiMeanPtr_();
377  }
378  else
379  {
380  return phiPtr_();
381  }
382 }
383 
384 
386 {
388  {
389  Info<< "Restoring field values to initial ones" << endl;
390  pInst() == pInitPtr_();
391  UInst() == UInitPtr_();
392  phiInst() == phiInitPtr_();
393  RASModelVariables_().restoreInitValues();
394  }
395 }
396 
397 
399 {
400  if (solverControl_.average())
401  {
402  Info<< "Resetting mean fields to zero" << endl;
403 
404  // Reset fields to zero
405  pMeanPtr_() == dimensionedScalar(pInst().dimensions(), Zero);
406  UMeanPtr_() == dimensionedVector(UInst().dimensions(), Zero);
407  phiMeanPtr_() == dimensionedScalar(phiInst().dimensions(), Zero);
408  RASModelVariables_().resetMeanFields();
410  // Reset averaging iteration index to 0
412  }
413 }
414 
415 
417 {
419  {
420  Info<< "Averaging fields" << endl;
421  label& iAverageIter = solverControl_.averageIter();
422  scalar avIter(iAverageIter);
423  scalar oneOverItP1 = 1./(avIter + 1);
424  scalar mult = avIter*oneOverItP1;
425  pMeanPtr_() == pMeanPtr_()*mult + pInst()*oneOverItP1;
426  UMeanPtr_() == UMeanPtr_()*mult + UInst()*oneOverItP1;
427  phiMeanPtr_() == phiMeanPtr_()*mult + phiInst()*oneOverItP1;
428  RASModelVariables_().computeMeanFields();
429  ++iAverageIter;
430  }
431 }
432 
433 
435 {
437  RASModelVariables_().correctBoundaryConditions(turbulence_());
438 }
439 
442 {
444 }
445 
448 {
449  return solverControl_.average();
450 }
451 
452 
454 {
455  incompressibleVars& incoVars = refCast<incompressibleVars>(vars);
456  // Copy source fields to the ones known by the object
457  swapAndRename(pPtr_, incoVars.pPtr_);
458  swapAndRename(UPtr_, incoVars.UPtr_);
459  swapAndRename(phiPtr_, incoVars.phiPtr_);
461  // Transfer turbulent fields. Copies fields since original fields are
462  // not owned by RASModelVariables but from the turbulence model
463  RASModelVariables_->transfer(incoVars.RASModelVariables()());
464 }
465 
466 
467 bool incompressibleVars::write() const
468 {
469  // Write dummy fields, for continuation only
471  {
472  if (RASModelVariables_().hasTMVar1())
473  {
474  createZeroFieldPtr<scalar>
475  (
476  mesh_,
477  RASModelVariables_().TMVar1BaseName(),
478  RASModelVariables_().TMVar1Inst().dimensions()
479  )().write();
480  }
481  if (RASModelVariables_().hasTMVar2())
482  {
483  createZeroFieldPtr<scalar>
484  (
485  mesh_,
486  RASModelVariables_().TMVar2BaseName(),
487  RASModelVariables_().TMVar2Inst().dimensions()
488  )().write();
489  }
490  if (RASModelVariables_().hasNut())
491  {
492  createZeroFieldPtr<scalar>
493  (
494  mesh_,
495  RASModelVariables_().nutBaseName(),
496  RASModelVariables_().nutRefInst().dimensions()
497  )().write();
498  }
499 
500  return true;
501  }
502 
503  return false;
504 }
505 
506 
507 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
508 
509 } // End namespace Foam
510 
511 // ************************************************************************* //
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:76
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:72
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.
Request registration (bool: true)
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