algebraic.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 "algebraic.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace XiModels
34 {
35  defineTypeNameAndDebug(algebraic, 0);
36  addToRunTimeSelectionTable(XiModel, algebraic, dictionary);
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const dictionary& XiProperties,
46  const psiuReactionThermo& thermo,
48  const volScalarField& Su,
49  const volScalarField& rho,
50  const volScalarField& b,
51  const surfaceScalarField& phi
52 )
53 :
54  XiModel(XiProperties, thermo, turbulence, Su, rho, b, phi),
55  XiShapeCoef(XiModelCoeffs_.lookup<scalar>("XiShapeCoef")),
56  XiEqModel_(XiEqModel::New(XiProperties, thermo, turbulence, Su)),
57  XiGModel_(XiGModel::New(XiProperties, thermo, turbulence, Su))
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
62 
64 {}
65 
66 
67 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
68 
70 {
71  return XiGModel_->Db();
72 }
73 
74 
76 {
77  volScalarField XiEqEta(XiEqModel_->XiEq());
78  volScalarField GEta(XiGModel_->G());
79 
80  volScalarField R(GEta*XiEqEta/(XiEqEta - 0.999));
81 
82  volScalarField XiEqStar(R/(R - GEta));
83 
84  Xi_ == 1.0 + (1.0 + (2*XiShapeCoef)*(0.5 - b_))*(XiEqStar - 1.0);
85 }
86 
87 
88 bool Foam::XiModels::algebraic::read(const dictionary& XiProperties)
89 {
90  XiModel::read(XiProperties);
91 
92  XiModelCoeffs_.lookup("XiShapeCoef") >> XiShapeCoef;
93 
94  return true;
95 }
96 
97 
98 // ************************************************************************* //
virtual tmp< volScalarField > Db() const
Return the flame diffusivity.
fluidReactionThermo & thermo
Definition: createFields.H:28
const volScalarField & b_
Definition: XiModel.H:121
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
volScalarField Xi_
Flame wrinkling field.
Definition: XiModel.H:125
Macros for easy insertion into run-time selection tables.
virtual bool read(const dictionary &XiProperties)
Update properties from given dictionary.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
algebraic(const dictionary &XiProperties, const psiuReactionThermo &thermo, const compressible::RASModel &turbulence, const volScalarField &Su, const volScalarField &rho, const volScalarField &b, const surfaceScalarField &phi)
Construct from components.
stressControl lookup("compactNormalStress") >> compactNormalStress
virtual bool read(const dictionary &XiProperties)=0
Update properties from given dictionary.
phi
Definition: correctPhi.H:3
dictionary XiModelCoeffs_
Definition: XiModel.H:115
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
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
#define R(A, B, C, D, E, F, K, M)
virtual ~algebraic()
Destructor.
const tmp< volScalarField::Internal > & Su
Definition: alphaSuSp.H:6
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
RASModel< momentumTransportModel > RASModel
Typedefs for turbulence, RAS and LES models for compressible flow based on the standard laminar trans...
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864
virtual void correct()
Correct the flame-wrinkling Xi.