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-2018 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 "pureMixture.H"
30 #include "multiComponentMixture.H"
31 #include "rhoThermo.H"
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
35 template<class Thermo, class OtherThermo>
36 template<class ThermoType>
39 (
40  const word& speciesName,
41  const multiComponentMixture<ThermoType>& globalThermo
42 ) const
43 {
44  return
45  globalThermo.getLocalThermo
46  (
47  globalThermo.species()
48  [
49  speciesName
50  ]
51  );
52 }
53 
54 
55 template<class Thermo, class OtherThermo>
56 template<class ThermoType>
59 (
60  const word& speciesName,
61  const pureMixture<ThermoType>& globalThermo
62 ) const
63 {
64  return globalThermo.cellMixture(0);
65 }
66 
67 
68 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
69 
70 template<class Thermo, class OtherThermo>
72 (
73  const dictionary& dict,
74  const phasePair& pair
75 )
76 :
77  interfaceCompositionModel(dict, pair),
78  thermo_
79  (
80  pair.phase1().mesh().lookupObject<Thermo>
81  (
82  IOobject::groupName(basicThermo::dictName, pair.phase1().name())
83  )
84  ),
85  otherThermo_
86  (
87  pair.phase2().mesh().lookupObject<OtherThermo>
88  (
89  IOobject::groupName(basicThermo::dictName, pair.phase2().name())
90  )
91  ),
92  Le_("Le", dimless, dict)
93 {}
94 
95 
96 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
97 
98 template<class Thermo, class OtherThermo>
101 {}
102 
103 
104 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
105 
106 template<class Thermo, class OtherThermo>
109 (
110  const word& speciesName,
111  const volScalarField& Tf
112 ) const
113 {
114  return
115  Yf(speciesName, Tf)
116  - thermo_.composition().Y()
117  [
118  thermo_.composition().species()[speciesName]
119  ];
120 }
121 
122 
123 template<class Thermo, class OtherThermo>
126 (
127  const word& speciesName
128 ) const
129 {
130  const typename Thermo::thermoType& localThermo =
132  (
133  speciesName,
134  thermo_
135  );
136 
137  const volScalarField& p(thermo_.p());
138 
139  const volScalarField& T(thermo_.T());
140 
141  tmp<volScalarField> tmpD
142  (
144  (
146  p.mesh(),
148  )
149  );
150 
151  volScalarField& D(tmpD.ref());
152 
153  forAll(p, celli)
154  {
155  D[celli] =
156  localThermo.alphah(p[celli], T[celli])
157  /localThermo.rho(p[celli], T[celli]);
158  }
159 
160  D /= Le_;
161 
162  return tmpD;
163 }
164 
165 
166 template<class Thermo, class OtherThermo>
169 (
170  const word& speciesName,
171  const volScalarField& Tf
172 ) const
173 {
174  const typename Thermo::thermoType& localThermo =
176  (
177  speciesName,
178  thermo_
179  );
180  const typename OtherThermo::thermoType& otherLocalThermo =
182  (
183  speciesName,
185  );
186 
187  const volScalarField& p(thermo_.p());
188  const volScalarField& otherP(otherThermo_.p());
189 
190  tmp<volScalarField> tmpL
191  (
193  (
195  p.mesh(),
197  )
198  );
199 
200  volScalarField& L(tmpL.ref());
201 
202  forAll(p, celli)
203  {
204  L[celli] =
205  localThermo.Ha(p[celli], Tf[celli])
206  - otherLocalThermo.Ha(otherP[celli], Tf[celli]);
207  }
208 
209  return tmpL;
210 }
211 
212 
213 template<class Thermo, class OtherThermo>
215 (
216  const volScalarField& K,
217  const volScalarField& Tf,
218  volScalarField& mDotL,
219  volScalarField& mDotLPrime
220 ) const
221 {
222  forAllConstIter(hashedWordList, this->speciesNames_, iter)
223  {
224  volScalarField rhoKDL
225  (
226  thermo_.rhoThermo::rho()
227  *K
228  *D(*iter)
229  *L(*iter, Tf)
230  );
231 
232  mDotL += rhoKDL*dY(*iter, Tf);
233  mDotLPrime += rhoKDL*YfPrime(*iter, Tf);
234  }
235 }
236 
237 
238 // ************************************************************************* //
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const OtherThermo & otherThermo_
Other Thermo.
virtual void addMDotL(const volScalarField &K, const volScalarField &Tf, volScalarField &mDotL, volScalarField &mDotLPrime) const
Add latent heat flow rate to total.
virtual tmp< volScalarField > dY(const word &speciesName, const volScalarField &Tf) const
Mass fraction difference between the interface and the field.
virtual tmp< volScalarField > D(const word &speciesName) const
Mass diffusivity.
virtual word name() const
Pair name.
~InterfaceCompositionModel()
Destructor.
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const
The interface species fraction.
virtual tmp< volScalarField > L(const word &speciesName, const volScalarField &Tf) const
Latent heat.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
const pureMixture< ThermoType >::thermoType & getLocalThermo(const word &speciesName, const pureMixture< ThermoType > &globalThermo) const
Get a reference to the local thermo for a pure mixture.
dynamicFvMesh & mesh
const phasePair & pair_
Phase pair.
static word groupName(Name name, const word &group)
const dimensionedScalar Le_
Lewis number.
const word dictName("particleTrackDict")
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
phaseModel & phase1
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const dimensionSet dimEnergy
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
const hashedWordList speciesNames_
Names of the transferring species.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
Definition: pureMixture.H:64
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual tmp< volScalarField > YfPrime(const word &speciesName, const volScalarField &Tf) const
The interface species fraction derivative w.r.t. temperature.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
phaseModel & phase2
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
InterfaceCompositionModel(const dictionary &dict, const phasePair &pair)
Construct from components.
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57