linearEquilibrium_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 
27 #include "laminarFlameSpeed.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace SuModels
35 {
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
43 
45 {
46  return SuModel::readCoeffs(dict);
47 
48  sigmaExt_.read(dict);
49 }
50 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
55 (
56  const dictionary& dict,
59 )
60 :
62  sigmaExt_("sigmaExt", dimless/dimTime, dict)
63 {
65 }
66 
67 
68 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
69 
71 {}
72 
73 
74 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
75 
77 {
78  const fvMesh& mesh(thermo_.mesh());
79 
80  const volVectorField& U(turbulence_.U());
83 
84  const volScalarField sigmas
85  (
86  ((n & n)*fvc::div(U) - (n & fvc::grad(U) & n))/Xi
87  + (
88  (n & n)*fvc::div(Su_*n)
89  - (n & fvc::grad(Su_*n) & n)
90  )*(Xi + scalar(1))/(2*Xi)
91  );
92 
93  Su_ == Su0_()()*max(scalar(1) - sigmas/sigmaExt_, scalar(0.01));
94 }
95 
96 
97 // ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
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
Equilibrium model for strained laminar flame speed with linear dependence on the strain-rate.
virtual void correct()
Correct the laminar flame speed.
virtual bool readCoeffs(const dictionary &dict)
Update coefficients from given dictionary.
linearEquilibrium(const dictionary &dict, const psiuMulticomponentThermo &thermo, const fluidThermoThermophysicalTransportModel &turbulence)
Construct from components.
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
void read(const dictionary &, const unitConversion &defaultUnits=NullObjectRef< unitConversion >())
Update the value of dimensioned<Type>
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
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::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
U
Definition: pEqn.H:72
addToRunTimeSelectionTable(SuModel, linearEquilibrium, dictionary)
defineTypeNameAndDebug(linearEquilibrium, 0)
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
Namespace for OpenFOAM.
const dimensionSet dimless
const dimensionSet dimTime
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dictionary dict
autoPtr< incompressible::momentumTransportModel > turbulence(incompressible::momentumTransportModel::New(U, phi, viscosity))
fluidMulticomponentThermo & thermo
Definition: createFields.H:31