Raoult.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 
26 #include "Raoult.H"
27 #include "phasePair.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace interfaceCompositionModels
35 {
36  defineTypeNameAndDebug(Raoult, 0);
37  addToRunTimeSelectionTable(interfaceCompositionModel, Raoult, dictionary);
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 (
46  const dictionary& dict,
47  const phasePair& pair
48 )
49 :
50  interfaceCompositionModel(dict, pair),
51  YNonVapour_
52  (
53  IOobject
54  (
55  IOobject::groupName("YNonVapour", pair.name()),
56  pair.phase1().mesh().time().timeName(),
57  pair.phase1().mesh()
58  ),
59  pair.phase1().mesh(),
61  ),
62  YNonVapourPrime_
63  (
64  IOobject
65  (
66  IOobject::groupName("YNonVapourPrime", pair.name()),
67  pair.phase1().mesh().time().timeName(),
68  pair.phase1().mesh()
69  ),
70  pair.phase1().mesh(),
72  )
73 {
74  forAllConstIter(hashedWordList, species(), iter)
75  {
76  speciesModels_.insert
77  (
78  *iter,
79  autoPtr<interfaceCompositionModel>
80  (
82  (
83  dict.subDict(*iter),
84  pair
85  )
86  )
87  );
88  }
89 }
90 
91 
92 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
93 
95 {}
96 
97 
98 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
99 
101 {
102  YNonVapour_ = scalar(1);
103 
104  forAllIter
105  (
106  HashTable<autoPtr<interfaceCompositionModel>>,
107  speciesModels_,
108  iter
109  )
110  {
111  iter()->update(Tf);
112 
113  YNonVapour_ -=
114  otherComposition().Y(iter.key())
115  *iter()->Yf(iter.key(), Tf);
116 
117  YNonVapourPrime_ -=
118  otherComposition().Y(iter.key())
119  *iter()->YfPrime(iter.key(), Tf);
120  }
121 }
122 
123 
125 (
126  const word& speciesName,
127  const volScalarField& Tf
128 ) const
129 {
130  if (species().contains(speciesName))
131  {
132  return
133  otherComposition().Y(speciesName)
134  *speciesModels_[speciesName]->Yf(speciesName, Tf);
135  }
136  else
137  {
138  return composition().Y(speciesName)*YNonVapour_;
139  }
140 }
141 
142 
145 (
146  const word& speciesName,
147  const volScalarField& Tf
148 ) const
149 {
150  if (species().contains(speciesName))
151  {
152  return
153  otherComposition().Y(speciesName)
154  *speciesModels_[speciesName]->YfPrime(speciesName, Tf);
155  }
156  else
157  {
158  return composition().Y(speciesName)*YNonVapourPrime_;
159  }
160 }
161 
162 
163 // ************************************************************************* //
const hashedWordList & species() const
Return the transferring species names.
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
static autoPtr< interfaceCompositionModel > New(const dictionary &dict, const phasePair &pair)
const basicSpecieMixture & otherComposition() const
Return the other composition.
Macros for easy insertion into run-time selection tables.
const phasePair & pair() const
Return the phase pair.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dynamicFvMesh & mesh
virtual void update(const volScalarField &Tf)
Update the composition.
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
word timeName
Definition: getTimeIndex.H:3
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const
The interface species fraction.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
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
Raoult(const dictionary &dict, const phasePair &pair)
Construct from components.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:53
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
virtual tmp< volScalarField > YfPrime(const word &speciesName, const volScalarField &Tf) const
The interface species fraction derivative w.r.t. temperature.
Namespace for OpenFOAM.