XiFluid.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) 2022-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 "XiFluid.H"
27 #include "localEulerDdtScheme.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace solvers
35 {
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 :
47  (
48  mesh,
50  ),
51 
53 
54  b_(thermo_.Y("b")),
55 
57 
58  Su
59  (
60  IOobject
61  (
62  "Su",
63  runTime.name(),
64  mesh,
65  IOobject::MUST_READ,
66  IOobject::AUTO_WRITE
67  ),
68  mesh
69  ),
70 
71  SuMin(0.01*Su.average()),
72  SuMax(4*Su.average()),
73 
74  Xi_
75  (
76  IOobject
77  (
78  "Xi",
79  runTime.name(),
80  mesh,
81  IOobject::MUST_READ,
82  IOobject::AUTO_WRITE
83  ),
84  mesh
85  ),
86 
87  St
88  (
89  IOobject
90  (
91  "St",
92  runTime.name(),
93  mesh,
94  IOobject::NO_READ,
95  IOobject::AUTO_WRITE
96  ),
97  Xi_*Su
98  ),
99 
100  combustionProperties
101  (
102  IOobject
103  (
104  "combustionProperties",
105  runTime.constant(),
106  mesh,
107  IOobject::MUST_READ_IF_MODIFIED,
108  IOobject::NO_WRITE
109  )
110  ),
111 
112  SuModel
113  (
114  combustionProperties.lookup("SuModel")
115  ),
116 
117  sigmaExt("sigmaExt", dimless/dimTime, combustionProperties),
118 
119  XiModel
120  (
121  combustionProperties.lookup("XiModel")
122  ),
123 
124  XiCoef("XiCoef", dimless, combustionProperties),
125 
126  XiShapeCoef("XiShapeCoef", dimless, combustionProperties),
127 
128  uPrimeCoef("uPrimeCoef", dimless, combustionProperties),
129 
130  ign(combustionProperties, runTime, mesh),
131 
133  (
134  momentumTransport(),
135  thermo_,
136  true
137  ),
138 
139  thermo(thermo_),
140  b(b_),
141  Xi(Xi_)
142 {
143  thermo.validate(type(), "ha", "ea");
144 
145  if (thermo_.containsSpecie("ft"))
146  {
147  fields.add(thermo_.Y("ft"));
148  }
149 
150  fields.add(b);
151  fields.add(thermo.he());
152  fields.add(thermo.heu());
153 }
154 
155 
156 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
157 
159 {}
160 
161 
162 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
163 
165 {
167 
168  if (pimple.predictTransport())
169  {
170  thermophysicalTransport.predict();
171  }
172 }
173 
174 
176 {
178 
179  if (pimple.correctTransport())
180  {
181  thermophysicalTransport.correct();
182  }
183 }
184 
185 
186 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
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
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
void validate(const string &app, const word &) const
Check that the thermodynamics package is consistent.
Definition: basicThermo.C:308
virtual const volScalarField & he() const =0
Enthalpy/Internal energy [J/kg].
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:57
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Abstract class for laminar flame speed.
Base-class for combustion fluid thermodynamic properties based on compressibility.
virtual const volScalarField & heu() const =0
Unburnt gas enthalpy [J/kg].
virtual PtrList< volScalarField > & Y()=0
Access the mass-fraction fields.
bool containsSpecie(const word &specieName) const
Does the mixture include this specie?
Abstract base class for run-time selectable region solvers.
Definition: solver.H:56
Solver module for compressible premixed/partially-premixed combustion with turbulence modelling.
Definition: XiFluid.H:97
const volScalarField & b
Reference to the combustion regress variable.
Definition: XiFluid.H:217
virtual void prePredictor()
Called at the start of the PIMPLE loop.
Definition: XiFluid.C:164
XiFluid(fvMesh &mesh)
Construct from region mesh.
Definition: XiFluid.C:44
psiuMulticomponentThermo & thermo_
Definition: XiFluid.H:103
virtual void postCorrector()
Correct the momentum and thermophysical transport modelling.
Definition: XiFluid.C:175
virtual ~XiFluid()
Destructor.
Definition: XiFluid.C:158
const psiuMulticomponentThermo & thermo
Reference to the fluid thermophysical properties.
Definition: XiFluid.H:213
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Set of fields used for the multivariate convection scheme.
Definition: XiFluid.H:112
Solver module for steady or transient turbulent flow of compressible isothermal fluids with optional ...
virtual void prePredictor()
Called at the start of the PIMPLE loop.
Definition: prePredictor.C:30
virtual void postCorrector()
Correct the momentum and thermophysical transport modelling.
pimpleControl pimple(mesh)
Info<< "Creating field dpdt\n"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.name(), mesh), mesh, dimensionedScalar(p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));Info<< "Creating the unstrained laminar flame speed\n"<< endl;autoPtr< laminarFlameSpeed > unstrainedLaminarFlameSpeed(laminarFlameSpeed::New(thermo))
volScalarField & b
Definition: createFields.H:25
Info<< "Creating thermophysical transport model\n"<< endl;turbulenceThermophysicalTransportModels::unityLewisEddyDiffusivity< RASThermophysicalTransportModel< ThermophysicalTransportModel< compressibleMomentumTransportModel, fluidThermo > >> thermophysicalTransport(turbulence(), thermo, true)
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
tmp< VolField< Type > > average(const SurfaceField< Type > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:46
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
addToRunTimeSelectionTable(solver, compressibleMultiphaseVoF, fvMesh)
defineTypeNameAndDebug(compressibleMultiphaseVoF, 0)
Namespace for OpenFOAM.
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:129
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
const dimensionSet dimTime
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
PtrList< volScalarField > & Y
fluidMulticomponentThermo & thermo
Definition: createFields.H:31