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 | 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 "XiEqModel.H"
27 
28 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32  defineTypeNameAndDebug(XiEqModel, 0);
33  defineRunTimeSelectionTable(XiEqModel, dictionary);
34 }
35 
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
39 Foam::XiEqModel::XiEqModel
40 (
41  const dictionary& XiEqProperties,
42  const psiuReactionThermo& thermo,
44  const volScalarField& Su
45 )
46 :
47  XiEqModelCoeffs_
48  (
49  XiEqProperties.subDict
50  (
51  word(XiEqProperties.lookup("XiEqModel")) + "Coeffs"
52  )
53  ),
54  thermo_(thermo),
55  turbulence_(turbulence),
56  Su_(Su)
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
61 
63 {}
64 
65 
66 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
67 
68 bool Foam::XiEqModel::read(const dictionary& XiEqProperties)
69 {
70  XiEqModelCoeffs_ = XiEqProperties.subDict(type() + "Coeffs");
71 
72  return true;
73 }
74 
75 
77 {
78  //***HGW It is not clear why B is written here
79  if (Su_.mesh().foundObject<volSymmTensorField>("B"))
80  {
81  const volSymmTensorField& B =
82  Su_.mesh().lookupObject<volSymmTensorField>("B");
83  B.write();
84  }
85 }
86 
87 
89 Foam::XiEqModel::calculateSchelkinEffect(const scalar uPrimeCoef) const
90 {
91  const fvMesh& mesh = Su_.mesh();
92 
93  const volVectorField& U = mesh.lookupObject<volVectorField>("U");
94  const volSymmTensorField& CT = mesh.lookupObject<volSymmTensorField>("CT");
95  const volScalarField& Nv = mesh.lookupObject<volScalarField>("Nv");
96  const volSymmTensorField& nsv =
97  mesh.lookupObject<volSymmTensorField>("nsv");
98 
99  tmp<volScalarField> tN
100  (
101  new volScalarField
102  (
103  IOobject
104  (
105  "tN",
106  mesh.time().timeName(),
107  mesh,
110  false
111  ),
112  mesh,
113  dimensionedScalar("zero", Nv.dimensions(), 0.0)
114  )
115  );
116  volScalarField& N = tN.ref();
117  N.primitiveFieldRef() = Nv.primitiveField()*pow(mesh.V(), 2.0/3.0);
118 
120  (
121  IOobject
122  (
123  "tns",
124  mesh.time().timeName(),
125  mesh,
128  ),
129  mesh,
131  (
132  "zero",
133  nsv.dimensions(),
134  Zero
135  )
136  );
137  ns.primitiveFieldRef() = nsv.primitiveField()*pow(mesh.V(), 2.0/3.0);
138 
139  const volVectorField Uhat
140  (
141  U/(mag(U) + dimensionedScalar("Usmall", U.dimensions(), 1e-4))
142  );
143 
144  const volScalarField nr(sqrt(max(N - (Uhat & ns & Uhat), scalar(1e-4))));
145 
146  const scalarField cellWidth(pow(mesh.V(), 1.0/3.0));
147 
148  const scalarField upLocal(uPrimeCoef*sqrt((U & CT & U)*cellWidth));
149 
150  const scalarField deltaUp(upLocal*(max(scalar(1.0), pow(nr, 0.5)) - 1.0));
151 
152  // Re use tN
153  N.primitiveFieldRef() = upLocal*(max(scalar(1.0), pow(nr, 0.5)) - 1.0);
154 
155  return tN;
156 }
157 
158 
159 // ************************************************************************* //
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
dictionary XiEqModelCoeffs_
Dictionary.
Definition: XiEqModel.H:65
const double e
Elementary charge.
Definition: doubleFloat.H:78
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:633
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual bool read(const dictionary &XiEqProperties)=0
Update properties from given dictionary.
tmp< volScalarField > calculateSchelkinEffect(const scalar) const
Return the sub-grid Schelkin effect.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
virtual ~XiEqModel()
Destructor.
void writeFields() const
Write fields.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
psiReactionThermo & thermo
Definition: createFields.H:31
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static const zero Zero
Definition: zero.H:91
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
Internal & ref()
Return a reference to the dimensioned internal field.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const Mesh & mesh() const
Return mesh.
virtual bool write() const
Write using setting from DB.
dimensioned< scalar > mag(const dimensioned< Type > &)
const volScalarField & Su_
Laminar burning velocity.
Definition: XiEqModel.H:74
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.