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-2024 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 "XiEqModel.H"
28 #include "XiProfile.H"
29 #include "XiGModel.H"
30 #include "fvmDiv.H"
31 #include "fvcLaplacian.H"
32 #include "fvModels.H"
33 #include "fvConstraints.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace XiModels
41 {
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
49 
51 {
53 
54  return true;
55 }
56 
57 
58 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
59 
61 (
62  const dictionary& dict,
65  const volScalarField& Su
66 )
67 :
69  XiEqModel_(XiEqModel::New(dict, thermo, turbulence, Su)),
70  XiProfile_(XiProfile::New(dict, b_)),
71  XiGModel_(XiGModel::New(dict, thermo, turbulence, Su))
72 {
75 }
76 
77 
78 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
79 
81 {}
82 
83 
84 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
85 
87 {
88  return XiEqModel_->Db();
89 }
90 
91 
93 {
94  const fvMesh& mesh(thermo_.mesh());
95 
98 
99  const volScalarField XiEqEta(XiEqModel_->XiEq());
100  const volScalarField GEta(XiGModel_->G());
101 
102  const volScalarField R(GEta*XiEqEta/(XiEqEta - 0.999));
103 
104  const volScalarField XiEqStar(R/(R - GEta));
105 
106  const volScalarField XiEq(1 + XiProfile_->profile()*(XiEqStar - 1));
107 
108  const volScalarField G(R*(XiEq - 1)/XiEq);
109 
110  const volScalarField& mgb = mesh.lookupObject<volScalarField>("mgb");
111  const surfaceScalarField& phiSt =
116 
117  const surfaceScalarField phiXi
118  (
119  "phiXi",
120  phiSt
121  + (
122  - fvc::interpolate(fvc::laplacian(Db, b_)/mgb)*nf
123  + fvc::interpolate(rho_)*fvc::interpolate(Su_*(1/Xi_ - Xi_))*nf
124  )
125  );
126 
127  const surfaceScalarField& phi = turbulence_.alphaRhoPhi();
128 
129  const volVectorField& U(turbulence_.U());
130 
131  const volVectorField Ut(U + Su_*Xi_*n);
132  const volScalarField sigmat((n & n)*fvc::div(Ut) - (n & fvc::grad(Ut) & n));
133 
134  const volScalarField sigmas
135  (
136  ((n & n)*fvc::div(U) - (n & fvc::grad(U) & n))/Xi_
137  + (
138  (n & n)*fvc::div(Su_*n)
139  - (n & fvc::grad(Su_*n) & n)
140  )*(Xi_ + scalar(1))/(2*Xi_)
141  );
142 
143  fvScalarMatrix XiEqn
144  (
145  fvm::ddt(rho_, Xi_)
146  + fvm::div(phi + phiXi, Xi_, "div(phiXi,Xi)")
147  - fvm::Sp(fvc::div(phiXi), Xi_)
148  ==
149  rho_*R
150  - fvm::Sp(rho_*(R - G), Xi_)
151  - fvm::Sp
152  (
153  rho_*max
154  (
155  sigmat - sigmas,
156  dimensionedScalar(sigmat.dimensions(), 0)
157  ),
158  Xi_
159  )
160  + fvModels.source(rho_, Xi_)
161  );
162 
163  XiEqn.relax();
164 
165  fvConstraints.constrain(XiEqn);
166 
167  XiEqn.solve();
168 
170 
171  // Limit range of Xi for realisability and stability
172  Xi_.max(1);
173  Xi_ = min(Xi_, 2*XiEq);
174 }
175 
176 
177 // ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
volScalarField Db("Db", rho *turbulence->nuEff())
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
writeOption writeOpt() const
Definition: IOobject.H:367
Templated abstract base class for thermophysical transport models.
Base class for equilibrium flame wrinkling XiEq models.
Definition: XiEqModel.H:52
Base class for flame wrinkling Xi generation rate models.
Definition: XiGModel.H:68
Base-class for all Xi models used by the b-Xi combustion model.
Definition: XiModel.H:117
volScalarField Xi_
Flame wrinkling field.
Definition: XiModel.H:139
virtual bool readCoeffs(const dictionary &dict)=0
Update coefficients from given dictionary.
Definition: XiModel.C:39
Transport model for flame wrinkling Xi.
Definition: transport.H:74
transport(const dictionary &dict, const psiuMulticomponentThermo &thermo, const fluidThermoThermophysicalTransportModel &turbulence, const volScalarField &Su)
Construct from components.
Definition: transport.C:61
virtual void correct()
Correct the flame-wrinkling Xi.
Definition: transport.C:92
virtual bool readCoeffs(const dictionary &dict)
Update coefficients from given dictionary.
Definition: transport.C:50
virtual tmp< volScalarField > Db() const
Return the flame diffusivity.
Definition: transport.C:86
virtual ~transport()
Destructor.
Definition: transport.C:80
Base class for flame wrinkling profiles.
Definition: XiProfile.H:50
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:103
Finite volume constraints.
Definition: fvConstraints.H:67
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:603
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
Finite volume models.
Definition: fvModels.H:65
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
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
Foam::fvConstraints & fvConstraints(Foam::fvConstraints::New(mesh))
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
Calculate the laplacian of the given field.
Calculate the matrix for the divergence of the given field and flux.
U
Definition: pEqn.H:72
defineTypeNameAndDebug(equilibrium, 0)
addToRunTimeSelectionTable(XiModel, equilibrium, dictionary)
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< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
tmp< DimensionedField< TypeR, GeoMesh, Field > > New(const tmp< DimensionedField< TypeR, GeoMesh, Field >> &tdf1, const word &name, const dimensionSet &dimensions)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict
autoPtr< incompressible::momentumTransportModel > turbulence(incompressible::momentumTransportModel::New(U, phi, viscosity))
fluidMulticomponentThermo & thermo
Definition: createFields.H:31