interfaceCompositionModel.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) 2015-2022 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 "phaseModel.H"
28 #include "phaseSystem.H"
29 #include "rhoReactionThermo.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  defineTypeNameAndDebug(interfaceCompositionModel, 0);
36  defineSidedInterfacialModelTypeNameAndDebug(interfaceCompositionModel, 0);
37  defineRunTimeSelectionTable(interfaceCompositionModel, dictionary);
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const dictionary& dict,
46  const phaseInterface& interface
47 )
48 :
49  interface_
50  (
51  interface.modelCast<interfaceCompositionModel, sidedPhaseInterface>()
52  ),
53  species_(dict.lookup("species")),
54  Le_("Le", dimless, dict),
55  thermo_(refCast<const rhoReactionThermo>(interface_.phase().thermo())),
56  otherThermo_(interface_.otherPhase().thermo())
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
61 
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67 
69 (
70  const word& speciesName,
71  const volScalarField& Tf
72 ) const
73 {
74  const label speciei = composition().species()[speciesName];
75 
76  return Yf(speciesName, Tf) - composition().Y()[speciei];
77 }
78 
79 
81 (
82  const word& speciesName,
83  const volScalarField& Tf
84 ) const
85 {
86  return YfPrime(speciesName, Tf);
87 }
88 
89 
91 (
92  const word& speciesName
93 ) const
94 {
95  const label speciei = composition().species()[speciesName];
96  const volScalarField& p(thermo_.p());
97  const volScalarField& T(thermo_.T());
98 
99  return volScalarField::New
100  (
101  IOobject::groupName("D" + speciesName, interface_.name()),
102  composition().kappa(speciei, p, T)
103  /composition().Cp(speciei, p, T)
104  /composition().rho(speciei, p, T)
105  /Le_
106  );
107 }
108 
109 
110 // ************************************************************************* //
dictionary dict
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
fluidReactionThermo & thermo
Definition: createFields.H:28
tmp< volScalarField > D(const word &speciesName) const
Mass diffusivity.
virtual scalar kappa(const label speciei, const scalar p, const scalar T) const =0
Thermal conductivity [W/m/K].
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:106
tmp< volScalarField > dYfPrime(const word &speciesName, const volScalarField &Tf) const
Mass fraction difference between the interface and the field.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
const dimensionSet dimless
virtual scalar rho(const label speciei, const scalar p, const scalar T) const =0
Density [kg/m^3].
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const
The interface species fraction.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
stressControl lookup("compactNormalStress") >> compactNormalStress
virtual scalar Cp(const label speciei, const scalar p, const scalar T) const =0
Heat capacity at constant pressure [J/kg/K].
interfaceCompositionModel(const dictionary &dict, const phaseInterface &interface)
Construct from a dictionary and an interface.
static word groupName(Name name, const word &group)
virtual ~interfaceCompositionModel()
Destructor.
#define defineSidedInterfacialModelTypeNameAndDebug(ModelType, DebugSwitch)
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const basicSpecieMixture & composition() const
Return the composition.
virtual tmp< volScalarField > YfPrime(const word &speciesName, const volScalarField &Tf) const
The interface species fraction derivative w.r.t. temperature.
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
tmp< volScalarField > dY(const word &speciesName, const volScalarField &Tf) const
Mass fraction difference between the interface and the field.
Namespace for OpenFOAM.
const speciesTable & species() const
Return the table of species.