algebraic.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-2015 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 "algebraic.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace XiModels
36 {
37  defineTypeNameAndDebug(algebraic, 0);
38  addToRunTimeSelectionTable(XiModel, algebraic, dictionary);
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
45 Foam::XiModels::algebraic::algebraic
46 (
47  const dictionary& XiProperties,
48  const psiuReactionThermo& thermo,
50  const volScalarField& Su,
51  const volScalarField& rho,
52  const volScalarField& b,
53  const surfaceScalarField& phi
54 )
55 :
56  XiModel(XiProperties, thermo, turbulence, Su, rho, b, phi),
57  XiShapeCoef(XiModelCoeffs_.get<scalar>("XiShapeCoef")),
58  XiEqModel_(XiEqModel::New(XiProperties, thermo, turbulence, Su)),
59  XiGModel_(XiGModel::New(XiProperties, thermo, turbulence, Su))
60 {}
61 
62 
63 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
64 
66 {}
67 
68 
69 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
70 
72 {
73  return XiGModel_->Db();
74 }
75 
76 
78 {
79  volScalarField XiEqEta(XiEqModel_->XiEq());
80  volScalarField GEta(XiGModel_->G());
81 
82  volScalarField R(GEta*XiEqEta/(XiEqEta - 0.999));
83 
84  volScalarField XiEqStar(R/(R - GEta));
85 
86  Xi_ == 1.0 + (1.0 + (2*XiShapeCoef)*(0.5 - b_))*(XiEqStar - 1.0);
87 }
88 
89 
90 bool Foam::XiModels::algebraic::read(const dictionary& XiProperties)
91 {
92  XiModel::read(XiProperties);
93 
94  XiModelCoeffs_.readEntry("XiShapeCoef", XiShapeCoef);
95 
96  return true;
97 }
98 
99 
100 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
virtual tmp< volScalarField > Db() const
Return the flame diffusivity.
zeroField Su
Definition: alphaSuSp.H:1
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Macros for easy insertion into run-time selection tables.
virtual bool read(const dictionary &XiProperties)
Update properties from given dictionary.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
psiReactionThermo & thermo
Definition: createFields.H:28
virtual bool read(const dictionary &XiProperties)=0
Update properties from given dictionary.
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
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)
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
#define R(A, B, C, D, E, F, K, M)
virtual ~algebraic()
Destructor.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
virtual void correct()
Correct the flame-wrinkling Xi.