transport_SuModel.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) 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_SuModel.H"
27 #include "laminarFlameSpeed.H"
28 #include "fvmDiv.H"
29 #include "fvcLaplacian.H"
30 #include "fvModels.H"
31 #include "fvConstraints.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace SuModels
39 {
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
47 
49 {
51 
52  sigmaExt_.read(dict);
53 
54  return true;
55 }
56 
57 
58 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
59 
61 (
62  const dictionary& dict,
65 )
66 :
68  sigmaExt_("sigmaExt", dimless/dimTime, dict),
69  SuMin_(0.01*Su_.average()),
70  SuMax_(4*Su_.average())
71 {
73 }
74 
75 
76 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
77 
79 {}
80 
81 
82 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
83 
85 {
86  Su_ = Su0_()();
87 }
88 
89 
91 {
92  const fvMesh& mesh(thermo_.mesh());
93 
96 
97  const volScalarField& rho = turbulence_.rho();
98  const volScalarField& b = thermo_.Y("b");
99  const volScalarField& mgb = mesh.lookupObject<volScalarField>("mgb");
100  const volScalarField& Xi = mesh.lookupObject<volScalarField>("Xi");
101  const surfaceScalarField& phiSt =
106 
107  const 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/Xi - Xi))*nf
114  )
115  );
116 
117  const surfaceScalarField& phi = turbulence_.alphaRhoPhi();
118  const volVectorField& U(turbulence_.U());
119 
120  const volScalarField sigmas
121  (
122  ((n & n)*fvc::div(U) - (n & fvc::grad(U) & n))/Xi
123  + (
124  (n & n)*fvc::div(Su_*n)
125  - (n & fvc::grad(Su_*n) & n)
126  )*(Xi + scalar(1))/(2*Xi)
127  );
128 
129  const volScalarField Su0(Su0_()());
130 
131  const volScalarField SuInf
132  (
133  Su0*max(scalar(1) - sigmas/sigmaExt_, scalar(0.01))
134  );
135 
136  const volScalarField Rc
137  (
138  (sigmas*SuInf*(Su0 - SuInf) + sqr(SuMin_)*sigmaExt_)
139  /(sqr(Su0 - SuInf) + sqr(SuMin_))
140  );
141 
142  fvScalarMatrix SuEqn
143  (
144  fvm::ddt(rho, Su_)
145  + fvm::div(phi + phiXi, Su_, "div(phiXi,Su)")
146  - fvm::Sp(fvc::div(phiXi), Su_)
147  ==
148  - fvm::SuSp(-rho*Rc*Su0/Su_, Su_)
149  - fvm::SuSp(rho*(sigmas + Rc), Su_)
150  + fvModels.source(rho, Su_)
151  );
152 
153  SuEqn.relax();
154 
155  fvConstraints.constrain(SuEqn);
156 
157  SuEqn.solve();
158 
160 
161  // Limit the maximum and minimum Su
162  Su_.min(SuMax_);
163  Su_.max(SuMin_);
164 }
165 
166 
167 // ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
volScalarField Db("Db", rho *turbulence->nuEff())
Generic GeometricField class.
writeOption writeOpt() const
Definition: IOobject.H:367
Base-class for all Su models used by the b-Xi combustion model.
Definition: SuModel.H:65
volScalarField Su_
Laminar flame speed.
Definition: SuModel.H:81
virtual bool readCoeffs(const dictionary &dict)=0
Update coefficients from given dictionary.
Definition: SuModel.C:39
Transport model for the strained laminar flame speed.
transport(const dictionary &dict, const psiuMulticomponentThermo &thermo, const fluidThermoThermophysicalTransportModel &turbulence)
Construct from components.
virtual void correct()
Correct the laminar flame speed.
virtual bool readCoeffs(const dictionary &dict)
Update coefficients from given dictionary.
virtual void reset()
Reset Su to the unburnt state.
virtual ~transport()
Destructor.
Simple unstrained model for the laminar flame speed.
Definition: unstrained.H:56
Templated abstract base class for thermophysical transport models.
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
void read(const dictionary &, const unitConversion &defaultUnits=NullObjectRef< unitConversion >())
Update the value of dimensioned<Type>
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.
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.
volScalarField & b
Definition: createFields.H:27
U
Definition: pEqn.H:72
addToRunTimeSelectionTable(SuModel, linearEquilibrium, dictionary)
defineTypeNameAndDebug(linearEquilibrium, 0)
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< 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 > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
const dimensionSet dimless
dimensioned< Type > average(const DimensionedField< Type, GeoMesh, PrimitiveField > &df)
const dimensionSet dimTime
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
dictionary dict
autoPtr< incompressible::momentumTransportModel > turbulence(incompressible::momentumTransportModel::New(U, phi, viscosity))
fluidMulticomponentThermo & thermo
Definition: createFields.H:31