basicXiSubG.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 "basicXiSubG.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace XiGModels
36 {
37  defineTypeNameAndDebug(basicSubGrid, 0);
38  addToRunTimeSelectionTable(XiGModel, basicSubGrid, dictionary);
39 };
40 };
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
45 Foam::XiGModels::basicSubGrid::basicSubGrid
46 (
47  const dictionary& XiGProperties,
48  const psiuReactionThermo& thermo,
50  const volScalarField& Su
51 )
52 :
53  XiGModel(XiGProperties, thermo, turbulence, Su),
54  k1(XiGModelCoeffs_.get<scalar>("k1")),
55  XiGModel_(XiGModel::New(XiGModelCoeffs_, thermo, turbulence, Su))
56 {}
57 
58 
59 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
60 
62 {}
63 
64 
65 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
66 
68 {
69  const objectRegistry& db = Su_.db();
70  const volVectorField& U = db.lookupObject<volVectorField>("U");
71  const volScalarField& Nv = db.lookupObject<volScalarField>("Nv");
72  const volScalarField& Lobs = db.lookupObject<volScalarField>("Lobs");
73 
74  tmp<volScalarField> tGtot = XiGModel_->G();
75  volScalarField& Gtot = tGtot.ref();
76 
77  const scalarField Cw = pow(Su_.mesh().V(), 2.0/3.0);
78  scalarField N(Nv.primitiveField()*Cw);
79 
80  forAll(N, celli)
81  {
82  if (N[celli] > 1e-3)
83  {
84  Gtot[celli] += k1*mag(U[celli])/Lobs[celli];
85  }
86  }
87 
88  return tGtot;
89 }
90 
91 
93 {
94  const objectRegistry& db = Su_.db();
95  const volScalarField& Xi = db.lookupObject<volScalarField>("Xi");
96  const volScalarField& rho = db.lookupObject<volScalarField>("rho");
97  const volScalarField& mgb = db.lookupObject<volScalarField>("mgb");
98  const volScalarField& Lobs = db.lookupObject<volScalarField>("Lobs");
99 
100  return XiGModel_->Db()
101  + rho*Su_*(Xi - 1.0)*mgb*(0.5*Lobs)*Lobs/(mgb*Lobs + 1.0);
102 }
103 
104 
105 bool Foam::XiGModels::basicSubGrid::read(const dictionary& XiGProperties)
106 {
107  XiGModel::read(XiGProperties);
108 
109  XiGModelCoeffs_.readEntry("k1", k1);
110 
111  return true;
112 }
113 
114 
115 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
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::tmp
A class for managing temporary objects.
Definition: PtrList.H:57
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::XiGModels::basicSubGrid::read
virtual bool read(const dictionary &XiGProperties)
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
Su
zeroField Su
Definition: alphaSuSp.H:1
Foam::XiGModels::basicSubGrid::Db
virtual tmp< volScalarField > Db() const
Foam::XiGModel::read
virtual bool read(const dictionary &XiGProperties)=0
rho
rho
Definition: readInitialConditions.H:88
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
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
basicXiSubG.H
Foam
Definition: atmBoundaryLayer.C:26
Foam::XiGModels::basicSubGrid::G
virtual tmp< volScalarField > G() const
Foam::XiGModels::basicSubGrid::~basicSubGrid
virtual ~basicSubGrid()
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::GeometricField::ref
Internal & ref(const bool updateAccessTime=true)
Definition: GeometricField.C:742
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Definition: createFields.H:11
Foam::PtrListOps::get
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
N
const Vector< label > N(dict.get< Vector< label >>("N"))
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)