basicXiSubXiEq.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 | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "basicXiSubXiEq.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace XiEqModels
34 {
35  defineTypeNameAndDebug(basicSubGrid, 0);
36  addToRunTimeSelectionTable(XiEqModel, basicSubGrid, dictionary);
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
43 Foam::XiEqModels::basicSubGrid::basicSubGrid
44 (
45  const dictionary& XiEqProperties,
46  const psiuReactionThermo& thermo,
48  const volScalarField& Su
49 )
50 :
51  XiEqModel(XiEqProperties, thermo, turbulence, Su),
52 
53  B_
54  (
55  IOobject
56  (
57  "B",
58  Su.mesh().facesInstance(),
59  Su.mesh(),
60  IOobject::MUST_READ,
61  IOobject::NO_WRITE
62  ),
63  Su.mesh()
64  ),
65 
66  XiEqModel_(XiEqModel::New(XiEqModelCoeffs_, thermo, turbulence, Su))
67 {}
68 
69 
70 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
71 
73 {}
74 
75 
76 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
77 
79 {
80  const fvMesh& mesh = Su_.mesh();
81  const volVectorField& U = mesh.lookupObject<volVectorField>("U");
82 
83  const volScalarField& Nv = mesh.lookupObject<volScalarField>("Nv");
84  const volSymmTensorField& nsv =
85  mesh.lookupObject<volSymmTensorField>("nsv");
86 
87  volScalarField magU(mag(U));
88  volVectorField Uhat
89  (
90  U/(mag(U) + dimensionedScalar("Usmall", U.dimensions(), 1e-4))
91  );
92 
93  const scalarField Cw = pow(mesh.V(), 2.0/3.0);
94 
96  (
97  IOobject
98  (
99  "N",
100  mesh.time().constant(),
101  mesh,
104  ),
105  mesh,
106  dimensionedScalar("zero", Nv.dimensions(), 0.0)
107  );
108  N.primitiveFieldRef() = Nv.primitiveField()*Cw;
109 
111  (
112  IOobject
113  (
114  "ns",
115  U.mesh().time().timeName(),
116  U.mesh(),
119  ),
120  U.mesh(),
122  (
123  "zero",
124  nsv.dimensions(),
125  Zero
126  )
127  );
128  ns.primitiveFieldRef() = nsv.primitiveField()*Cw;
129 
130  volScalarField n(max(N - (Uhat & ns & Uhat), scalar(1e-4)));
131  volScalarField b((Uhat & B_ & Uhat)/sqrt(n));
132  volScalarField up(sqrt((2.0/3.0)*turbulence_.k()));
133 
134  volScalarField XiSubEq
135  (
136  scalar(1)
137  + max(2.2*sqrt(b), min(0.34*magU/up*sqrt(b), scalar(1.6)))
138  * min(n, scalar(1))
139  );
140 
141  return (XiSubEq*XiEqModel_->XiEq());
142 }
143 
144 
145 bool Foam::XiEqModels::basicSubGrid::read(const dictionary& XiEqProperties)
146 {
147  XiEqModel::read(XiEqProperties);
148 
149  return XiEqModel_->read(XiEqModelCoeffs_);
150 }
151 
152 
153 // ************************************************************************* //
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
const double e
Elementary charge.
Definition: doubleFloat.H:78
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual bool read(const dictionary &XiEqProperties)=0
Update properties from given dictionary.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
Macros for easy insertion into run-time selection tables.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
psiReactionThermo & thermo
Definition: createFields.H:31
dynamicFvMesh & mesh
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static const zero Zero
Definition: zero.H:91
virtual tmp< volScalarField > XiEq() const
Return the flame-wrinking XiEq.
virtual ~basicSubGrid()
Destructor.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
virtual bool read(const dictionary &XiEqProperties)
Update properties from given dictionary.
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label N
Definition: createFields.H:22
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
A class for managing temporary objects.
Definition: PtrList.H:54
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:44
Namespace for OpenFOAM.