transport.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-2023 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 "transport.H"
27 #include "surfaceInterpolate.H"
28 #include "fvmDdt.H"
29 #include "fvcLaplacian.H"
30 #include "fvmDiv.H"
31 #include "fvmSup.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace XiModels
39 {
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
49 (
50  const dictionary& XiProperties,
53  const volScalarField& Su,
54  const volScalarField& rho,
55  const volScalarField& b,
56  const surfaceScalarField& phi
57 )
58 :
59  XiModel(XiProperties, thermo, turbulence, Su, rho, b, phi),
60  XiShapeCoef(XiModelCoeffs_.lookup<scalar>("XiShapeCoef")),
61  XiEqModel_(XiEqModel::New(XiProperties, thermo, turbulence, Su)),
62  XiGModel_(XiGModel::New(XiProperties, thermo, turbulence, Su))
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
67 
69 {}
70 
71 
72 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
73 
75 {
76  return XiGModel_->Db();
77 }
78 
79 
81 (
83 )
84 {
85  volScalarField XiEqEta(XiEqModel_->XiEq());
86  volScalarField GEta(XiGModel_->G());
87 
88  volScalarField R(GEta*XiEqEta/(XiEqEta - 0.999));
89 
90  volScalarField XiEqStar(R/(R - GEta));
91 
92  volScalarField XiEq
93  (
94  1.0 + (1.0 + (2*XiShapeCoef)*(0.5 - b_))*(XiEqStar - 1.0)
95  );
96 
97  volScalarField G(R*(XiEq - 1.0)/XiEq);
98 
99  const objectRegistry& db = b_.db();
100  const volScalarField& betav = db.lookupObject<volScalarField>("betav");
101  const volScalarField& mgb = db.lookupObject<volScalarField>("mgb");
102  const surfaceScalarField& phiSt =
103  db.lookupObject<surfaceScalarField>("phiSt");
104  const volScalarField& Db = db.lookupObject<volScalarField>("Db");
105  const surfaceScalarField& nf = db.lookupObject<surfaceScalarField>("nf");
106 
107  surfaceScalarField phiXi
108  (
109  "phiXi",
110  phiSt
111  + (
112  - fvc::interpolate(fvc::laplacian(Db, b_)/mgb)*nf
113  + fvc::interpolate(rho_)*fvc::interpolate(Su_*(1.0/Xi_ - Xi_))*nf
114  )
115  );
116 
117  solve
118  (
119  betav*fvm::ddt(rho_, Xi_)
120  + mvConvection.fvmDiv(phi_, Xi_)
121  + fvm::div(phiXi, Xi_)
122  - fvm::Sp(fvc::div(phiXi), Xi_)
123  ==
124  betav*rho_*R
125  - fvm::Sp(betav*rho_*(R - G), Xi_)
126  );
127 
128  // Correct boundedness of Xi
129  // ~~~~~~~~~~~~~~~~~~~~~~~~~
130  Xi_.max(1.0);
131  Xi_ = min(Xi_, 2.0*XiEq);
132 }
133 
134 
136 {
137  XiModel::read(XiProperties);
138 
139  XiModelCoeffs_.lookup("XiShapeCoef") >> XiShapeCoef;
140 
141  return true;
142 }
143 
144 
145 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
volScalarField Db("Db", rho *turbulence->nuEff())
tmp< fv::convectionScheme< scalar > > mvConvection(fv::convectionScheme< scalar >::New(mesh, fields, phi, mesh.schemes().div("div(phi,ft_b_ha_hau)")))
Generic GeometricField class.
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:312
Base-class for all XiEq models used by the b-XiEq combustion model. The available models are : basicX...
Definition: XiEqModel.H:58
Base-class for all Xi generation models used by the b-Xi combustion model. See Technical Report SH/RE...
Definition: XiGModel.H:55
Base-class for all Xi models used by the b-Xi combustion model. See Technical Report SH/RE/01R for de...
Definition: XiModel.H:109
virtual bool read(const dictionary &XiProperties)=0
Update properties from given dictionary.
Definition: XiModel.C:86
Simple transport model for Xi based on Gulders correlation with a linear correction function to give ...
Definition: transport.H:58
transport(const dictionary &XiProperties, const psiuMulticomponentThermo &thermo, const compressible::RASModel &turbulence, const volScalarField &Su, const volScalarField &rho, const volScalarField &b, const surfaceScalarField &phi)
Construct from components.
Definition: transport.C:49
virtual void correct()
Correct the flame-wrinkling Xi.
Definition: transport.H:110
virtual bool read(const dictionary &XiProperties)
Update properties from given dictionary.
Definition: transport.C:135
virtual tmp< volScalarField > Db() const
Return the flame diffusivity.
Definition: transport.C:74
virtual ~transport()
Destructor.
Definition: transport.C:68
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Abstract base class for convection schemes.
Registry of regIOobjects.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
Base-class for combustion fluid thermodynamic properties based on compressibility.
A class for managing temporary objects.
Definition: tmp.H:55
Calculate the laplacian of the given field.
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for implicit and explicit sources.
volScalarField & b
Definition: createFields.H:27
addToRunTimeSelectionTable(XiModel, algebraic, dictionary)
defineTypeNameAndDebug(algebraic, 0)
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
RASModel< momentumTransportModel > RASModel
Typedefs for turbulence, RAS and LES models for compressible flow based on the standard laminar trans...
const dimensionedScalar G
Newtonian constant of gravitation.
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvcLaplacian.C:45
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
autoPtr< incompressible::momentumTransportModel > turbulence(incompressible::momentumTransportModel::New(U, phi, viscosity))
fluidMulticomponentThermo & thermo
Definition: createFields.H:31