Loading...
Searching...
No Matches
XiEqModel.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) 2011-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "XiEqModel.H"
29
30// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31
32namespace Foam
33{
36}
37
38
39// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40
41Foam::XiEqModel::XiEqModel
42(
43 const dictionary& XiEqProperties,
44 const psiuReactionThermo& thermo,
45 const compressible::RASModel& turbulence,
46 const volScalarField& Su
47)
48:
49 XiEqModelCoeffs_
50 (
51 XiEqProperties.subDict
52 (
53 XiEqProperties.get<word>("XiEqModel") + "Coeffs"
54 )
55 ),
56 thermo_(thermo),
57 turbulence_(turbulence),
58 Su_(Su)
59{}
60
61
62// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
63
65{}
66
67
68// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
69
70bool Foam::XiEqModel::read(const dictionary& XiEqProperties)
71{
72 XiEqModelCoeffs_ = XiEqProperties.optionalSubDict(type() + "Coeffs");
73
74 return true;
75}
76
77
79{
80 //***HGW It is not clear why B is written here
81 const auto* B = Su_.mesh().cfindObject<volSymmTensorField>("B");
82 if (B)
83 {
84 B->write();
85 }
86}
87
88
89Foam::tmp<Foam::volScalarField>
90Foam::XiEqModel::calculateSchelkinEffect(const scalar uPrimeCoef) const
91{
92 const fvMesh& mesh = Su_.mesh();
93
97 const volSymmTensorField& nsv =
99
100 auto tN = volScalarField::New
101 (
102 "tN",
104 mesh,
105 dimensionedScalar(Nv.dimensions(), Zero)
106 );
107 auto& N = tN.ref();
108
109 N.primitiveFieldRef() = Nv.primitiveField()*pow(mesh.V(), 2.0/3.0);
110
111 auto tns = volSymmTensorField::New
112 (
113 "tns",
115 mesh,
116 dimensionedSymmTensor(nsv.dimensions(), Zero)
117 );
118 auto& ns = tns.ref();
119
120 ns.primitiveFieldRef() = nsv.primitiveField()*pow(mesh.V(), 2.0/3.0);
121
122 const volVectorField Uhat
123 (
124 U/(mag(U) + dimensionedScalar("Usmall", U.dimensions(), 1e-4))
125 );
126
127 const volScalarField nr(sqrt(max(N - (Uhat & ns & Uhat), scalar(1e-4))));
128
129 const scalarField cellWidth(pow(mesh.V(), 1.0/3.0));
130
131 const scalarField upLocal(uPrimeCoef*sqrt((U & CT & U)*cellWidth));
132
133 const scalarField deltaUp(upLocal*(max(scalar(1), pow(nr, 0.5)) - 1.0));
134
135 // Re use tN
136 N.primitiveFieldRef() = upLocal*(max(scalar(1), pow(nr, 0.5)) - 1.0);
137
138 return tN;
139}
140
141
142// ************************************************************************* //
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
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).
Base-class for all XiEq models used by the b-XiEq combustion model. The available models are : basicX...
Definition XiEqModel.H:56
virtual ~XiEqModel()
Destructor.
virtual bool read(const dictionary &XiEqProperties)=0
Update properties from given dictionary.
tmp< volScalarField > calculateSchelkinEffect(const scalar) const
Return the sub-grid Schelkin effect.
void writeFields() const
Write fields.
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
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
dynamicFvMesh & mesh
compressible::turbulenceModel & turbulence
zeroField Su
Definition alphaSuSp.H:1
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
volScalarField & e
const Vector< label > N(dict.get< Vector< label > >("N"))