multicomponentFluid.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 "multicomponentFluid.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  Y_(thermo_.Y()),
55 
56  reaction(combustionModel::New(thermo_, momentumTransport())),
57 
59  (
61  (
62  momentumTransport(),
63  thermo_
64  )
65  ),
66 
67  thermo(thermo_),
68  Y(Y_)
69 {
70  thermo.validate(type(), "h", "e");
71 
72  forAll(Y, i)
73  {
74  fields.add(Y[i]);
75  }
76  fields.add(thermo.he());
77 }
78 
79 
80 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
81 
83 {}
84 
85 
86 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
87 
89 {
91 
92  if (pimple.predictTransport())
93  {
94  thermophysicalTransport->predict();
95  }
96 }
97 
98 
100 {
102 
103  if (pimple.correctTransport())
104  {
105  thermophysicalTransport->correct();
106  }
107 }
108 
109 
110 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Templated abstract base class for thermophysical transport models.
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 combustion models.
Base-class for multi-component fluid thermodynamic properties.
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:57
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Reaction base-class holding the specie names and coefficients.
Definition: reaction.H:57
Abstract base class for run-time selectable region solvers.
Definition: solver.H:56
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.
Solver module for steady or transient turbulent flow of compressible multicomponent fluids with optio...
const PtrList< volScalarField > & Y
Reference to the composition.
virtual void prePredictor()
Called at the start of the PIMPLE loop.
virtual void postCorrector()
Correct the momentum and thermophysical transport modelling.
const fluidMulticomponentThermo & thermo
Reference to the fluid thermophysical properties.
multicomponentFluid(fvMesh &mesh)
Construct from region mesh.
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
pimpleControl pimple(mesh)
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)
addToRunTimeSelectionTable(solver, compressibleMultiphaseVoF, fvMesh)
defineTypeNameAndDebug(compressibleMultiphaseVoF, 0)
Namespace for OpenFOAM.
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:129
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