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-2019 PCOpt/NTUA
9 Copyright (C) 2013-2019 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
61(
62 fvMesh& mesh,
63 const word& managerType,
64 const dictionary& dict
65)
66:
68 solverControl_(SIMPLEControl::New(mesh, managerType, *this)),
69 incoVars_(allocateVars())
70{
72 (
77 );
78}
79
80
81// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
82
84{
85 const Time& time = mesh_.time();
86 Info<< "Time = " << time.timeName() << "\n" << endl;
87
88 // Grab references
90 incoVars_.turbulence();
91 turbulence->correct();
92
93 solverControl_().write();
94
95 // Average fields if necessary
96 incoVars_.computeMeanFields();
97
99}
100
101
103{
104 // Iterate
105 if (active_)
106 {
107 // Reset initial and mean fields before solving
108 while (solverControl_().loop())
109 {
110 solveIter();
111 }
112 }
113}
114
115
117{
118 return solverControl_().loop();
119}
120
121
123{
124 os.writeEntry("averageIter", solverControl_().averageIter());
125
126 return true;
127}
128
129
130// ************************************************************************* //
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
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
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.
Definition: SIMPLEControl.H:56
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Ostream & printExecutionTime(OSstream &os) const
Print the elapsed ExecutionTime (cpu-time), ClockTime.
Definition: TimeIO.C:618
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Base class for primal incompressible solvers.
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:82
virtual const dictionary & dict() const
Return the solver dictionary.
Definition: solver.C:113
fvMesh & mesh_
Reference to the mesh database.
Definition: solver.H:60
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
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,...
compressible::turbulenceModel & turbulence
Namespace for OpenFOAM.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
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:34
dictionary dict