XiEqModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 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 
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.optionalSubDict(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  (
102  (
103  "tN",
104  mesh,
105  dimensionedScalar(Nv.dimensions(), 0)
106  )
107  );
108  volScalarField& N = tN.ref();
109  N.primitiveFieldRef() = Nv.primitiveField()*pow(mesh.V(), 2.0/3.0);
110 
112  (
113  IOobject
114  (
115  "tns",
116  mesh.time().timeName(),
117  mesh,
120  ),
121  mesh,
123  (
124  "zero",
125  nsv.dimensions(),
126  Zero
127  )
128  );
129  ns.primitiveFieldRef() = nsv.primitiveField()*pow(mesh.V(), 2.0/3.0);
130 
131  const volVectorField Uhat
132  (
133  U/(mag(U) + dimensionedScalar(U.dimensions(), 1e-4))
134  );
135 
136  const volScalarField nr(sqrt(max(N - (Uhat & ns & Uhat), scalar(1e-4))));
137 
138  const scalarField cellWidth(pow(mesh.V(), 1.0/3.0));
139 
140  const scalarField upLocal(uPrimeCoef*sqrt((U & CT & U)*cellWidth));
141 
142  const scalarField deltaUp(upLocal*(max(scalar(1), pow(nr, 0.5)) - 1.0));
143 
144  // Re use tN
145  N.primitiveFieldRef() = upLocal*(max(scalar(1), pow(nr, 0.5)) - 1.0);
146 
147  return tN;
148 }
149 
150 
151 // ************************************************************************* //
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:62
dictionary XiEqModelCoeffs_
Dictionary.
Definition: XiEqModel.H:65
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
fluidReactionThermo & thermo
Definition: createFields.H:28
tmp< volScalarField > calculateSchelkinEffect(const scalar) const
Return the sub-grid Schelkin effect.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
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:59
fvMesh & mesh
virtual ~XiEqModel()
Destructor.
XiEqModel(const dictionary &XiEqProperties, const psiuReactionThermo &thermo, const compressible::RASModel &turbulence, const volScalarField &Su)
Construct from components.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
stressControl lookup("compactNormalStress") >> compactNormalStress
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:1060
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static const zero Zero
Definition: zero.H:97
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::momentumTransportModel > turbulence(compressible::momentumTransportModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
const Mesh & mesh() const
Return mesh.
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
Internal & ref()
Return a reference to the dimensioned internal field.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dimensioned< scalar > mag(const dimensioned< Type > &)
const volScalarField & Su_
Laminar burning velocity.
Definition: XiEqModel.H:74
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
const tmp< volScalarField::Internal > & Su
Definition: alphaSuSp.H:6
virtual bool write(const bool write=true) const
Write using setting from DB.
A class for managing temporary objects.
Definition: PtrList.H:53
void writeFields() const
Write fields.
RASModel< momentumTransportModel > RASModel
Typedefs for turbulence, RAS and LES models for compressible flow based on the standard laminar trans...
Namespace for OpenFOAM.