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-2020 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 "phasePair.H"
29 #include "phaseSystem.H"
30 #include "rhoReactionThermo.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(interfaceCompositionModel, 0);
37  defineRunTimeSelectionTable(interfaceCompositionModel, dictionary);
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const dictionary& dict,
46  const phasePair& pair
47 )
48 :
49  pair_(pair),
50  species_(dict.lookup("species")),
51  Le_("Le", dimless, dict),
52  phase_
53  (
54  // !!! This is a hack, to let the models know which side of the
55  // interface it applies to. In general, we are going to need improve
56  // the phase-pair system to take into account model sidedness and other
57  // types of asymmetry from just dispersed/continuous.
58  pair.phase1().name() == IOobject::group(dict.dictName())
59  ? pair.phase1()
60  : pair.phase2()
61  ),
62  otherPhase_(pair.otherPhase(phase_)),
63  thermo_
64  (
65  phase_.mesh().lookupObject<rhoReactionThermo>
66  (
67  IOobject::groupName(basicThermo::dictName, phase_.name())
68  )
69  ),
70  otherThermo_
71  (
72  otherPhase_.mesh().lookupObject<rhoThermo>
73  (
74  IOobject::groupName(basicThermo::dictName, otherPhase_.name())
75  )
76  )
77 {}
78 
79 
80 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
81 
83 {}
84 
85 
86 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
87 
89 (
90  const word& speciesName,
91  const volScalarField& Tf
92 ) const
93 {
94  const label speciei = composition().species()[speciesName];
95 
96  return Yf(speciesName, Tf) - composition().Y()[speciei];
97 }
98 
99 
101 (
102  const word& speciesName
103 ) const
104 {
105  const label speciei = composition().species()[speciesName];
106  const volScalarField& p(thermo_.p());
107  const volScalarField& T(thermo_.T());
108 
109  return volScalarField::New
110  (
111  IOobject::groupName("D", pair_.name()),
112  composition().alphah(speciei, p, T)
113  /composition().rho(speciei, p, T)
114  /Le_
115  );
116 }
117 
118 
120 (
121  const word& speciesName,
122  const volScalarField& Tf
123 ) const
124 {
125  const label speciei = composition().species()[speciesName];
126  const volScalarField& p(thermo_.p());
127  volScalarField Ha(composition().Ha(speciei, p, Tf));
128 
129  const volScalarField& otherP(otherThermo_.p());
130  tmp<volScalarField> otherHa(nullptr);
131  if (otherHasComposition())
132  {
133  const label otherSpeciei = otherComposition().species()[speciesName];
134  otherHa = otherComposition().Ha(otherSpeciei, otherP, Tf);
135  }
136  else
137  {
138  otherHa = otherThermo_.ha(otherP, Tf);
139  }
140 
141  return
143  (
144  IOobject::groupName("L", pair_.name()),
145  otherHa - Ha
146  );
147 }
148 
149 
151 (
152  const volScalarField& K,
153  const volScalarField& Tf,
154  volScalarField& dmdtL,
155  volScalarField& dmdtLPrime
156 ) const
157 {
158  forAllConstIter(hashedWordList, species_, iter)
159  {
160  const volScalarField rhoKDL(thermo_.rho()*K*D(*iter)*L(*iter, Tf));
161 
162  dmdtL += rhoKDL*dY(*iter, Tf);
163  dmdtLPrime += rhoKDL*YfPrime(*iter, Tf);
164  }
165 }
166 
167 
168 // ************************************************************************* //
const char *const group
Group name for atomic constants.
dictionary dict
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
interfaceCompositionModel(const dictionary &dict, const phasePair &pair)
Construct from a dictionary and a phase pair.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
virtual scalar Ha(const label speciei, const scalar p, const scalar T) const =0
Absolute enthalpy [J/kg].
virtual tmp< volScalarField > L(const word &speciesName, const volScalarField &Tf) const
Latent heat.
const basicSpecieMixture & otherComposition() const
Return the other composition.
virtual tmp< volScalarField > D(const word &speciesName) const
Mass diffusivity.
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.
virtual tmp< volScalarField > dY(const word &speciesName, const volScalarField &Tf) const
Mass fraction difference between the interface and the field.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
static word groupName(Name name, const word &group)
virtual ~interfaceCompositionModel()
Destructor.
const word dictName("particleTrackDict")
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
bool otherHasComposition() const
Return whether the other side has a multi-specie composition.
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
virtual scalar alphah(const label speciei, const scalar p, const scalar T) const =0
Thermal conductivity [W/m/K].
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
virtual void addDmdtL(const volScalarField &K, const volScalarField &Tf, volScalarField &dmdtL, volScalarField &dmdtLPrime) const
Add latent heat flow rate to total.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const basicSpecieMixture & composition() const
Return the composition.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
virtual tmp< volScalarField > YfPrime(const word &speciesName, const volScalarField &Tf) const
The interface species fraction derivative w.r.t. temperature.
scalar Ha(const scalar p, const scalar T) const
Definition: EtoHthermo.H:20
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
Namespace for OpenFOAM.
const speciesTable & species() const
Return the table of species.