RASTurbulenceModel.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-2023 PCOpt/NTUA
9  Copyright (C) 2013-2023 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 
30 #include "RASTurbulenceModel.H"
31 #include "findRefCell.H"
32 #include "Time.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(RASTurbulenceModel, 0);
41  (
42  incompressiblePrimalSolver,
43  RASTurbulenceModel,
45  );
46 }
47 
48 
49 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
50 
52 {
54  return getIncoVars();
55 }
56 
57 
58 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
59 
60 Foam::RASTurbulenceModel::RASTurbulenceModel
61 (
62  fvMesh& mesh,
63  const word& managerType,
64  const dictionary& dict,
65  const word& solverName
66 )
67 :
69  (
70  mesh,
71  managerType,
72  dict,
73  solverName
74  ),
75  solverControl_(SIMPLEControl::New(mesh, managerType, *this)),
76  incoVars_(allocateVars())
77 {
79  (
80  incoVars_.pInst(),
81  solverControl_().dict(),
82  solverControl_().pRefCell(),
83  solverControl_().pRefValue()
84  );
85 }
86 
87 
88 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
89 
91 {
92  const Time& time = mesh_.time();
93  Info<< "Time = " << time.timeName() << "\n" << endl;
94 
95  // Grab references
97  incoVars_.turbulence();
98  turbulence->correct();
99 
100  solverControl_().write();
101 
102  // Average fields if necessary
103  incoVars_.computeMeanFields();
104 
105  time.printExecutionTime(Info);
106 }
107 
108 
110 {
111  // Iterate
112  if (active_)
113  {
114  // Reset initial and mean fields before solving
115  while (solverControl_().loop())
116  {
117  solveIter();
118  }
119  }
120 }
121 
124 {
125  return solverControl_().loop();
126 }
127 
128 
130 {
131  os.writeEntry("averageIter", solverControl_().averageIter());
132 
133  return true;
134 }
135 
136 
137 // ************************************************************************* //
const volScalarField & pInst() const
Return const reference to pressure.
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
autoPtr< variablesSet > vars_
Base variableSet pointer.
Definition: solver.H:95
incompressibleVars & incoVars_
Reference to incompressibleVars.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Find the reference cell nearest (in index) to the given cell but which is not on a cyclic...
Macros for easy insertion into run-time selection tables.
Base class for solution control classes.
dynamicFvMesh & mesh
virtual void solve()
Main control loop.
A class for handling words, derived from Foam::string.
Definition: word.H:63
const Time & time() const noexcept
Return time registry.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
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
OBJstream os(runTime.globalPath()/outputName)
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
defineTypeNameAndDebug(combustionModel, 0)
Ostream & printExecutionTime(OSstream &os) const
Print the elapsed ExecutionTime (cpu-time), ClockTime.
Definition: TimeIO.C:604
virtual bool loop()
Looper (advances iters, time step)
virtual void solveIter()
Execute one iteration of the solution algorithm.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
incompressibleVars & allocateVars()
Protected Member Functions.
SIMPLE control class to supply convergence information/checks for the SIMPLE loop.
Definition: SIMPLEControl.H:46
messageStream Info
Information stream (stdout output on master, null elsewhere)
const incompressibleVars & getIncoVars() const
Access to the incompressible variables set.
fvMesh & mesh_
Reference to the mesh database.
Definition: solver.H:56
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Base class for primal incompressible solvers.
virtual bool writeData(Ostream &os) const
Read average iteration.
bool setRefCell(const volScalarField &field, const volScalarField &fieldRef, const dictionary &dict, label &refCelli, scalar &refValue, const bool forceReference=false)
If the field fieldRef needs referencing find the reference cell nearest.
Definition: findRefCell.C:27
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
autoPtr< SIMPLEControl > solverControl_
Solver control.