SpalartAllmaras.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-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 "SpalartAllmaras.H"
31 #include "wallDist.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace incompressible
39 {
40 namespace RASVariables
41 {
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
51 (
52  const fvMesh& mesh,
53  const solverControl& SolverControl
54 )
55 :
56  RASModelVariables(mesh, SolverControl)
57 {
58  TMVar1BaseName_ = "nuTilda";
59 
62 
63  // The wall dist name can vary depending on how wallDist was
64  // constructed. Grab the field directly from wallDist
65 
66  distPtr_.ref
67  (
68  const_cast<volScalarField&>(wallDist::New(mesh_).y())
69  );
70 
73 }
74 
75 
76 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
77 
79 (
81 ) const
82 {
83  auto tnutJacobian = tmp<volScalarField>::New
84  (
85  IOobject
86  (
87  "nutJacobianVar1",
88  mesh_.time().timeName(),
89  mesh_,
92  ),
93  mesh_,
95  );
96  auto& nutJacobian = tnutJacobian.ref();
97 
98  const volScalarField& nu = laminarTransport.nu();
99  const volScalarField& nuTilda = TMVar1();
100 
101  volScalarField chi(nuTilda/nu);
102  volScalarField chi3(pow3(chi));
103 
104  const scalar Cv13 = pow3(7.1);
105  volScalarField fv1(chi3/(chi3 + Cv13));
106  volScalarField fv1ByChi2Sqr(sqr(chi/(chi3 + Cv13)));
107  volScalarField Cdfv1(3.0*Cv13*fv1ByChi2Sqr);
108 
109  nutJacobian = Cdfv1*chi + fv1;
110 
111  return tnutJacobian;
112 }
113 
114 
115 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
116 
117 } // End namespace RASVariables
118 } // End namespace incompressible
119 } // End namespace Foam
120 
121 // ************************************************************************* //
Type & lookupObjectRef(const word &name, const bool recursive=false) const
Lookup and return non-const reference to the object of the given Type. Fatal if not found or the wron...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
static const wallDist & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
Macros for easy insertion into run-time selection tables.
virtual tmp< volScalarField > nutJacobianVar1(const singlePhaseTransportModel &laminarTransport) const
return nut Jacobian wrt the TM vars
scalar y
singlePhaseTransportModel laminarTransport(U, phi)
SpalartAllmaras(const fvMesh &mesh, const solverControl &SolverControl)
Construct from components.
dynamicFvMesh & mesh
Base class for solver control classes.
Definition: solverControl.H:45
const volScalarField & TMVar1() const
Return references to turbulence fields.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
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
dimensionedScalar pow3(const dimensionedScalar &ds)
addToRunTimeSelectionTable(RASModelVariables, kEpsilon, dictionary)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Nothing to be read.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
A simple single-phase transport model based on viscosityModel.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
volScalarField & nu
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127