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 -------------------------------------------------------------------------------
10 License
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 
32 namespace Foam
33 {
34  defineTypeNameAndDebug(XiEqModel, 0);
35  defineRunTimeSelectionTable(XiEqModel, dictionary);
36 }
37 
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
41 Foam::XiEqModel::XiEqModel
42 (
43  const dictionary& XiEqProperties,
44  const psiuReactionThermo& thermo,
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 
70 bool 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  if (Su_.mesh().foundObject<volSymmTensorField>("B"))
82  {
83  const volSymmTensorField& B =
84  Su_.mesh().lookupObject<volSymmTensorField>("B");
85  B.write();
86  }
87 }
88 
89 
91 Foam::XiEqModel::calculateSchelkinEffect(const scalar uPrimeCoef) const
92 {
93  const fvMesh& mesh = Su_.mesh();
94 
95  const volVectorField& U = mesh.lookupObject<volVectorField>("U");
96  const volSymmTensorField& CT = mesh.lookupObject<volSymmTensorField>("CT");
97  const volScalarField& Nv = mesh.lookupObject<volScalarField>("Nv");
98  const volSymmTensorField& nsv =
99  mesh.lookupObject<volSymmTensorField>("nsv");
100 
101  tmp<volScalarField> tN
102  (
103  new volScalarField
104  (
105  IOobject
106  (
107  "tN",
108  mesh.time().timeName(),
109  mesh,
113  ),
114  mesh,
115  dimensionedScalar(Nv.dimensions(), Zero)
116  )
117  );
118  volScalarField& N = tN.ref();
119  N.primitiveFieldRef() = Nv.primitiveField()*pow(mesh.V(), 2.0/3.0);
120 
122  (
123  IOobject
124  (
125  "tns",
126  mesh.time().timeName(),
127  mesh,
130  ),
131  mesh,
132  dimensionedSymmTensor(nsv.dimensions(), Zero)
133  );
134  ns.primitiveFieldRef() = nsv.primitiveField()*pow(mesh.V(), 2.0/3.0);
135 
136  const volVectorField Uhat
137  (
138  U/(mag(U) + dimensionedScalar("Usmall", U.dimensions(), 1e-4))
139  );
140 
141  const volScalarField nr(sqrt(max(N - (Uhat & ns & Uhat), scalar(1e-4))));
142 
143  const scalarField cellWidth(pow(mesh.V(), 1.0/3.0));
144 
145  const scalarField upLocal(uPrimeCoef*sqrt((U & CT & U)*cellWidth));
146 
147  const scalarField deltaUp(upLocal*(max(scalar(1), pow(nr, 0.5)) - 1.0));
148 
149  // Re use tN
150  N.primitiveFieldRef() = upLocal*(max(scalar(1), pow(nr, 0.5)) - 1.0);
151 
152  return tN;
153 }
154 
155 
156 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:85
zeroField Su
Definition: alphaSuSp.H:1
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
tmp< volScalarField > calculateSchelkinEffect(const scalar) const
Return the sub-grid Schelkin effect.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual bool read(const dictionary &XiEqProperties)=0
Update properties from given dictionary.
Ignore writing from objectRegistry::writeObject()
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
virtual ~XiEqModel()
Destructor.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
psiReactionThermo & thermo
Definition: createFields.H:28
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
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("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
const Vector< label > N(dict.get< Vector< label >>("N"))
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
U
Definition: pEqn.H:72
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Nothing to be read.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
void writeFields() const
Write fields.
Do not request registration (bool: false)
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133