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  composition(thermo_.composition()),
55 
56  b_(composition.Y("b")),
57 
59 
60  Su
61  (
62  IOobject
63  (
64  "Su",
65  runTime.name(),
66  mesh,
67  IOobject::MUST_READ,
68  IOobject::AUTO_WRITE
69  ),
70  mesh
71  ),
72 
73  SuMin(0.01*Su.average()),
74  SuMax(4*Su.average()),
75 
76  Xi_
77  (
78  IOobject
79  (
80  "Xi",
81  runTime.name(),
82  mesh,
83  IOobject::MUST_READ,
84  IOobject::AUTO_WRITE
85  ),
86  mesh
87  ),
88 
89  St
90  (
91  IOobject
92  (
93  "St",
94  runTime.name(),
95  mesh,
96  IOobject::NO_READ,
97  IOobject::AUTO_WRITE
98  ),
99  Xi_*Su
100  ),
101 
102  combustionProperties
103  (
104  IOobject
105  (
106  "combustionProperties",
107  runTime.constant(),
108  mesh,
109  IOobject::MUST_READ_IF_MODIFIED,
110  IOobject::NO_WRITE
111  )
112  ),
113 
114  SuModel
115  (
116  combustionProperties.lookup("SuModel")
117  ),
118 
119  sigmaExt
120  (
121  combustionProperties.lookup("sigmaExt")
122  ),
123 
124  XiModel
125  (
126  combustionProperties.lookup("XiModel")
127  ),
128 
129  XiCoef
130  (
131  combustionProperties.lookup("XiCoef")
132  ),
133 
134  XiShapeCoef
135  (
136  combustionProperties.lookup("XiShapeCoef")
137  ),
138 
139  uPrimeCoef
140  (
141  combustionProperties.lookup("uPrimeCoef")
142  ),
143 
144  ign(combustionProperties, runTime, mesh),
145 
147  (
148  momentumTransport(),
149  thermo_,
150  true
151  ),
152 
153  thermo(thermo_),
154  b(b_),
155  Xi(Xi_)
156 {
157  thermo.validate(type(), "ha", "ea");
158 
159  if (composition.contains("ft"))
160  {
161  fields.add(composition.Y("ft"));
162  }
163 
164  fields.add(b);
165  fields.add(thermo.he());
166  fields.add(thermo.heu());
167 }
168 
169 
170 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
171 
173 {}
174 
175 
176 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
177 
179 {
181 
182  if (pimple.predictTransport())
183  {
184  thermophysicalTransport.predict();
185  }
186 }
187 
188 
190 {
192 
193  if (pimple.correctTransport())
194  {
195  thermophysicalTransport.correct();
196  }
197 }
198 
199 
200 // ************************************************************************* //
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
volScalarField & Y(const word &specieName)
Return the mass-fraction field for a specie given by name.
bool contains(const word &specieName) const
Does the mixture include this specie?
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
void validate(const string &app, const word &) const
Check that the thermodynamics package is consistent.
Definition: basicThermo.C:331
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:57
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Abstract class for laminar flame speed.
Base-class for combustion fluid thermodynamic properties based on compressibility.
virtual volScalarField & heu()=0
Unburnt gas enthalpy [J/kg].
Abstract base class for run-time selectable region solvers.
Definition: solver.H:55
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:221
virtual void prePredictor()
Called at the start of the PIMPLE loop.
Definition: XiFluid.C:178
XiFluid(fvMesh &mesh)
Construct from region mesh.
Definition: XiFluid.C:44
virtual void postCorrector()
Correct the momentum and thermophysical transport modelling.
Definition: XiFluid.C:189
virtual ~XiFluid()
Destructor.
Definition: XiFluid.C:172
basicCombustionMixture & composition
Reference to the combustion mixture.
Definition: XiFluid.H:109
const psiuMulticomponentThermo & thermo
Reference to the fluid thermophysical properties.
Definition: XiFluid.H:217
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Set of fields used for the multivariate convection scheme.
Definition: XiFluid.H:116
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:27
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:111
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
basicSpecieMixture & composition
PtrList< volScalarField > & Y
fluidMulticomponentThermo & thermo
Definition: createFields.H:31