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-2015 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"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace XiEqModels
35 {
36  defineTypeNameAndDebug(basicSubGrid, 0);
37  addToRunTimeSelectionTable(XiEqModel, basicSubGrid, dictionary);
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
44 Foam::XiEqModels::basicSubGrid::basicSubGrid
45 (
46  const dictionary& XiEqProperties,
47  const psiuReactionThermo& thermo,
49  const volScalarField& Su
50 )
51 :
52  XiEqModel(XiEqProperties, thermo, turbulence, Su),
53 
54  B_
55  (
56  IOobject
57  (
58  "B",
59  Su.mesh().facesInstance(),
60  Su.mesh(),
61  IOobject::MUST_READ,
62  IOobject::NO_WRITE
63  ),
64  Su.mesh()
65  ),
66 
67  XiEqModel_(XiEqModel::New(XiEqModelCoeffs_, thermo, turbulence, Su))
68 {}
69 
70 
71 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
72 
74 {}
75 
76 
77 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
78 
80 {
81  const fvMesh& mesh = Su_.mesh();
82  const volVectorField& U = mesh.lookupObject<volVectorField>("U");
83 
84  const volScalarField& Nv = mesh.lookupObject<volScalarField>("Nv");
85  const volSymmTensorField& nsv =
86  mesh.lookupObject<volSymmTensorField>("nsv");
87 
88  volScalarField magU(mag(U));
89  volVectorField Uhat
90  (
91  U/(mag(U) + dimensionedScalar("Usmall", U.dimensions(), 1e-4))
92  );
93 
94  const scalarField Cw = pow(mesh.V(), 2.0/3.0);
95 
96  tmp<volScalarField> tN
97  (
98  new volScalarField
99  (
100  IOobject
101  (
102  "tN",
103  mesh.time().constant(),
104  mesh,
107  ),
108  mesh,
109  dimensionedScalar("zero", Nv.dimensions(), 0.0),
110  zeroGradientFvPatchVectorField::typeName
111  )
112  );
113 
114  volScalarField& N = tN();
115 
116  N.internalField() = Nv.internalField()*Cw;
117 
118  tmp<volSymmTensorField> tns
119  (
121  (
122  IOobject
123  (
124  "tns",
125  U.mesh().time().timeName(),
126  U.mesh(),
129  ),
130  U.mesh(),
132  (
133  "zero",
134  nsv.dimensions(),
135  pTraits<symmTensor>::zero
136  ),
137  zeroGradientFvPatchSymmTensorField::typeName
138  )
139  );
140 
141  volSymmTensorField& ns = tns();
142 
143  ns.internalField() = nsv.internalField()*Cw;
144 
145  volScalarField n(max(N - (Uhat & ns & Uhat), scalar(1e-4)));
146 
147  volScalarField b((Uhat & B_ & Uhat)/sqrt(n));
148 
149  volScalarField up(sqrt((2.0/3.0)*turbulence_.k()));
150 
151  volScalarField XiSubEq
152  (
153  scalar(1)
154  + max(2.2*sqrt(b), min(0.34*magU/up*sqrt(b), scalar(1.6)))
155  * min(n, scalar(1))
156  );
157 
158  return (XiSubEq*XiEqModel_->XiEq());
159 }
160 
161 
162 bool Foam::XiEqModels::basicSubGrid::read(const dictionary& XiEqProperties)
163 {
164  XiEqModel::read(XiEqProperties);
165 
166  return XiEqModel_->read(XiEqModelCoeffs_);
167 }
168 
169 
170 // ************************************************************************* //
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:28
dimensioned< scalar > mag(const dimensioned< Type > &)
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
dynamicFvMesh & mesh
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
virtual ~basicSubGrid()
Destructor.
Namespace for OpenFOAM.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:44
label n
const double e
Elementary charge.
Definition: doubleFloat.H:78
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Macros for easy insertion into run-time selection tables.
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
virtual tmp< volScalarField > XiEq() const
Return the flame-wrinking XiEq.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
virtual bool read(const dictionary &XiEqProperties)
Update properties from given dictionary.
virtual bool read(const dictionary &XiEqProperties)=0
Update properties from given dictionary.
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
psiReactionThermo & thermo
Definition: createFields.H:32
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:118
defineTypeNameAndDebug(combustionModel, 0)