Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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
30#include "RASTurbulenceModel.H"
31#include "findRefCell.H"
32#include "Time.H"
34
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
41 (
45 );
46}
47
48
49// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
50
52{
54 return getIncoVars();
55}
56
57
58// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
59
60Foam::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 (
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())
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// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const scalar pRefValue
const label pRefCell
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Solves for a RAS turbulence model, with constant U and phi values.
autoPtr< SIMPLEControl > solverControl_
Solver control.
virtual bool writeData(Ostream &os) const
Read average iteration.
incompressibleVars & allocateVars()
Protected Member Functions.
incompressibleVars & incoVars_
Reference to incompressibleVars.
virtual bool loop()
Looper (advances iters, time step).
virtual void solveIter()
Execute one iteration of the solution algorithm.
virtual void solve()
Main control loop.
SIMPLE control class to supply convergence information/checks for the SIMPLE loop.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
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
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Base class for primal incompressible solvers.
static autoPtr< incompressiblePrimalSolver > New(fvMesh &mesh, const word &managerType, const dictionary &dict, const word &solverName)
Return a reference to the selected incompressible primal solver.
const incompressibleVars & getIncoVars() const
Access to the incompressible variables set.
Base class for solution control classes.
const volScalarField & pInst() const
Return const reference to pressure.
const Time & time() const noexcept
Return time registry.
autoPtr< variablesSet > vars_
Base variableSet pointer.
Definition solver.H:95
const word & managerType() const
Return the manager type.
Definition solverI.H:72
bool active_
Solve equations?
Definition solver.H:76
const dictionary & dict() const
Return the solver dictionary.
Definition solverI.H:54
const fvMesh & mesh() const
Return the solver mesh.
Definition solverI.H:24
const word & solverName() const
Return the solver name.
Definition solverI.H:30
fvMesh & mesh_
Reference to the mesh database.
Definition solver.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
Find the reference cell nearest (in index) to the given cell but which is not on a cyclic,...
Namespace for OpenFOAM.
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.
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
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
dictionary dict