Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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 "SpalartAllmaras.H"
31#include "wallDist.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
40namespace 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
60 TMVar1Ptr_.ref(mesh_.lookupObjectRef<volScalarField>(TMVar1BaseName_));
61 nutPtr_.ref(mesh_.lookupObjectRef<volScalarField>(nutBaseName_));
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(wallDist::New(mesh_).y().constCast());
67
70}
71
72
73// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74
76(
78) const
79{
80 auto tnutJacobian = volScalarField::New
81 (
82 "nutJacobianVar1",
84 mesh_,
86 );
87 auto& nutJacobian = tnutJacobian.ref();
88
90 const volScalarField& nu = tnu();
91 const volScalarField& nuTilda = TMVar1();
92
93 volScalarField chi(nuTilda/nu);
94 volScalarField chi3(pow3(chi));
95
96 const scalar Cv13 = pow3(7.1);
97 volScalarField fv1(chi3/(chi3 + Cv13));
98 volScalarField fv1ByChi2Sqr(sqr(chi/(chi3 + Cv13)));
99 volScalarField Cdfv1(3.0*Cv13*fv1ByChi2Sqr);
100
101 nutJacobian = Cdfv1*chi + fv1;
102
103 return tnutJacobian;
104}
105
106
107// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
108
109} // End namespace RASVariables
110} // End namespace incompressible
111} // End namespace Foam
112
113// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
static FOAM_NO_DANGLING_REFERENCE const wallDist & New(const fvMesh &mesh, Args &&... args)
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
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
RASModelVariables(const fvMesh &mesh, const solverControl &SolverControl)
Construct from components.
const volScalarField & TMVar1() const
Return references to turbulence fields.
SpalartAllmaras(const fvMesh &mesh, const solverControl &SolverControl)
Construct from components.
virtual tmp< volScalarField > nutJacobianVar1(const singlePhaseTransportModel &laminarTransport) const
return nut Jacobian wrt the TM vars
A simple single-phase transport model based on viscosityModel.
Base class for solver control classes.
A class for managing temporary objects.
Definition tmp.H:75
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
volScalarField & nu
singlePhaseTransportModel laminarTransport(U, phi)