InterfaceCompositionModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2015-2018 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
29 #include "phaseModel.H"
30 #include "phasePair.H"
31 #include "pureMixture.H"
32 #include "multiComponentMixture.H"
33 #include "rhoThermo.H"
34 
35 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36 
37 template<class Thermo, class OtherThermo>
38 template<class ThermoType>
41 (
42  const word& speciesName,
43  const multiComponentMixture<ThermoType>& globalThermo
44 ) const
45 {
46  return
47  globalThermo.getLocalThermo
48  (
49  globalThermo.species()
50  [
51  speciesName
52  ]
53  );
54 }
55 
56 
57 template<class Thermo, class OtherThermo>
58 template<class ThermoType>
61 (
62  const word& speciesName,
63  const pureMixture<ThermoType>& globalThermo
64 ) const
65 {
66  return globalThermo.cellMixture(0);
67 }
68 
69 
70 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
71 
72 template<class Thermo, class OtherThermo>
74 (
75  const dictionary& dict,
76  const phasePair& pair
77 )
78 :
79  interfaceCompositionModel(dict, pair),
80  thermo_
81  (
82  pair.phase1().mesh().lookupObject<Thermo>
83  (
84  IOobject::groupName(basicThermo::dictName, pair.phase1().name())
85  )
86  ),
87  otherThermo_
88  (
89  pair.phase2().mesh().lookupObject<OtherThermo>
90  (
91  IOobject::groupName(basicThermo::dictName, pair.phase2().name())
92  )
93  ),
94  Le_("Le", dimless, dict)
95 {}
96 
97 
98 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
99 
100 template<class Thermo, class OtherThermo>
103 (
104  const word& speciesName,
105  const volScalarField& Tf
106 ) const
107 {
108  return
109  Yf(speciesName, Tf)
110  - thermo_.composition().Y()
111  [
112  thermo_.composition().species()[speciesName]
113  ];
114 }
115 
116 
117 template<class Thermo, class OtherThermo>
120 (
121  const word& speciesName
122 ) const
123 {
124  const typename Thermo::thermoType& localThermo =
125  getLocalThermo
126  (
127  speciesName,
128  thermo_
129  );
130 
131  const volScalarField& p(thermo_.p());
132 
133  const volScalarField& T(thermo_.T());
134 
135  tmp<volScalarField> tmpD
136  (
138  (
139  IOobject::groupName("D", pair_.name()),
140  p.mesh(),
142  )
143  );
144 
145  volScalarField& D = tmpD.ref();
146 
147  forAll(p, celli)
148  {
149  D[celli] =
150  localThermo.alphah(p[celli], T[celli])
151  /localThermo.rho(p[celli], T[celli]);
152  }
153 
154  D /= Le_;
155 
156  return tmpD;
157 }
158 
159 
160 template<class Thermo, class OtherThermo>
163 (
164  const word& speciesName,
165  const volScalarField& Tf
166 ) const
167 {
168  const typename Thermo::thermoType& localThermo =
169  getLocalThermo
170  (
171  speciesName,
172  thermo_
173  );
174  const typename OtherThermo::thermoType& otherLocalThermo =
175  getLocalThermo
176  (
177  speciesName,
178  otherThermo_
179  );
180 
181  const volScalarField& p(thermo_.p());
182  const volScalarField& otherP(otherThermo_.p());
183 
184  tmp<volScalarField> tmpL
185  (
187  (
188  IOobject::groupName("L", pair_.name()),
189  p.mesh(),
191  )
192  );
193 
194  volScalarField& L = tmpL.ref();
195 
196  forAll(p, celli)
197  {
198  L[celli] =
199  localThermo.Ha(p[celli], Tf[celli])
200  - otherLocalThermo.Ha(otherP[celli], Tf[celli]);
201  }
202 
203  return tmpL;
204 }
205 
206 
207 template<class Thermo, class OtherThermo>
209 (
210  const volScalarField& K,
211  const volScalarField& Tf,
212  volScalarField& mDotL,
213  volScalarField& mDotLPrime
214 ) const
215 {
216  for (const word& speciesName : this->speciesNames_)
217  {
218  volScalarField rhoKDL
219  (
220  thermo_.rhoThermo::rho()
221  *K
222  *D(speciesName)
223  *L(speciesName, Tf)
224  );
225 
226  mDotL += rhoKDL*dY(speciesName, Tf);
227  mDotLPrime += rhoKDL*YfPrime(speciesName, Tf);
228  }
229 }
230 
231 
232 // ************************************************************************* //
L
const vector L(dict.get< vector >("L"))
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::InterfaceCompositionModel::L
virtual tmp< volScalarField > L(const word &speciesName, const volScalarField &Tf) const
Latent heat (to - from)(thermo - otherThermo)
Definition: InterfaceCompositionModel.C:264
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::InterfaceCompositionModel::getLocalThermo
const pureMixture< ThermoType >::thermoType & getLocalThermo(const word &speciesName, const pureMixture< ThermoType > &globalThermo) const
Get a reference to the local thermo for a pure mixture.
dictName
const word dictName("faMeshDefinition")
Foam::pureMixture::thermoType
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
Definition: pureMixture.H:66
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::GeometricField< scalar, fvPatchField, volMesh >::New
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &ds, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field from name, mesh, dimensions and patch type.
Definition: GeometricFieldNew.C:34
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:58
rhoThermo.H
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::InterfaceCompositionModel::InterfaceCompositionModel
InterfaceCompositionModel(const dictionary &dict, const phasePair &pair)
Construct from components.
Definition: InterfaceCompositionModel.C:181
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dict
dictionary dict
Definition: searchingEngine.H:14
phase2
phaseModel & phase2
Definition: setRegionFluidFields.H:6
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
InterfaceCompositionModel.H
Foam::InterfaceCompositionModel::addMDotL
virtual void addMDotL(const volScalarField &K, const volScalarField &Tf, volScalarField &mDotL, volScalarField &mDotLPrime) const
Add latent heat flow rate to total.
Definition: InterfaceCompositionModel.C:209
multiComponentMixture.H
D
const dimensionedScalar & D
Definition: solveBulkSurfactant.H:4
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::multiComponentMixture::thermoType
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
Definition: multiComponentMixture.H:89
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
phase1
phaseModel & phase1
Definition: setRegionFluidFields.H:5
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::InterfaceCompositionModel::D
virtual tmp< volScalarField > D(const word &speciesName) const
Mass diffusivity of the local thermo.
Definition: InterfaceCompositionModel.C:218
pureMixture.H
Foam::InterfaceCompositionModel::dY
virtual tmp< volScalarField > dY(const word &speciesName, const volScalarField &Tf) const
Mass fraction difference between the interface and the field.
Definition: InterfaceCompositionModel.C:306
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189