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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 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 "basicXiSubXiEq.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace XiEqModels
36 {
37  defineTypeNameAndDebug(basicSubGrid, 0);
38  addToRunTimeSelectionTable(XiEqModel, basicSubGrid, dictionary);
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
45 Foam::XiEqModels::basicSubGrid::basicSubGrid
46 (
47  const dictionary& XiEqProperties,
48  const psiuReactionThermo& thermo,
50  const volScalarField& Su
51 )
52 :
53  XiEqModel(XiEqProperties, thermo, turbulence, Su),
54 
55  B_
56  (
57  IOobject
58  (
59  "B",
60  Su.mesh().facesInstance(),
61  Su.mesh(),
62  IOobject::MUST_READ,
63  IOobject::NO_WRITE
64  ),
65  Su.mesh()
66  ),
67 
68  XiEqModel_(XiEqModel::New(XiEqModelCoeffs_, thermo, turbulence, Su))
69 {}
70 
71 
72 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
73 
75 {}
76 
77 
78 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
79 
81 {
82  const fvMesh& mesh = Su_.mesh();
83  const volVectorField& U = mesh.lookupObject<volVectorField>("U");
84 
85  const volScalarField& Nv = mesh.lookupObject<volScalarField>("Nv");
86  const volSymmTensorField& nsv =
87  mesh.lookupObject<volSymmTensorField>("nsv");
88 
89  volScalarField magU(mag(U));
90  volVectorField Uhat
91  (
92  U/(mag(U) + dimensionedScalar("Usmall", U.dimensions(), 1e-4))
93  );
94 
95  const scalarField Cw = pow(mesh.V(), 2.0/3.0);
96 
98  (
99  IOobject
100  (
101  "N",
102  mesh.time().constant(),
103  mesh,
106  ),
107  mesh,
108  dimensionedScalar(Nv.dimensions(), Zero)
109  );
110  N.primitiveFieldRef() = Nv.primitiveField()*Cw;
111 
113  (
114  IOobject
115  (
116  "ns",
117  U.mesh().time().timeName(),
118  U.mesh(),
121  ),
122  U.mesh(),
123  dimensionedSymmTensor(nsv.dimensions(), Zero)
124  );
125  ns.primitiveFieldRef() = nsv.primitiveField()*Cw;
126 
127  volScalarField n(max(N - (Uhat & ns & Uhat), scalar(1e-4)));
128  volScalarField b((Uhat & B_ & Uhat)/sqrt(n));
129  volScalarField up(sqrt((2.0/3.0)*turbulence_.k()));
130 
131  volScalarField XiSubEq
132  (
133  scalar(1)
134  + max(2.2*sqrt(b), min(0.34*magU/up*sqrt(b), scalar(1.6)))
135  * min(n, scalar(1))
136  );
137 
138  return (XiSubEq*XiEqModel_->XiEq());
139 }
140 
141 
142 bool Foam::XiEqModels::basicSubGrid::read(const dictionary& XiEqProperties)
143 {
144  XiEqModel::read(XiEqProperties);
145 
146  return XiEqModel_->read(XiEqModelCoeffs_);
147 }
148 
149 
150 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:191
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:46
turbulence
Info<< "Reading field U\n"<< 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\n"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
Foam::dimensionedSymmTensor
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
Definition: dimensionedSymmTensor.H:46
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:57
Foam::Zero
static constexpr const zero Zero
Definition: zero.H:131
basicXiSubXiEq.H
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::XiEqModels::basicSubGrid::~basicSubGrid
virtual ~basicSubGrid()
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Definition: hashSets.C:26
Su
zeroField Su
Definition: alphaSuSp.H:1
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::XiEqModel::read
virtual bool read(const dictionary &XiEqProperties)=0
Foam::constant::physicoChemical::b
const dimensionedScalar b
Definition: createFields.H:27
Foam::volSymmTensorField
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:61
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:36
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:53
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:68
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Definition: hashSets.C:40
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Definition: atmBoundaryLayer.C:26
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Definition: DimensionedFieldReuseFunctions.H:100
U
U
Definition: pEqn.H:72
Foam::compressible::RASModel
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
Definition: turbulentFluidThermoModel.H:62
Foam::XiEqModels::basicSubGrid::XiEq
virtual tmp< volScalarField > XiEq() const
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:137
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Definition: createFields.H:11
N
const Vector< label > N(dict.get< Vector< label >>("N"))
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:184
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::XiEqModels::basicSubGrid::read
virtual bool read(const dictionary &XiEqProperties)