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-2023 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"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace interfaceCompositionModels
34 {
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const dictionary& dict,
46  const phaseInterface& interface
47 )
48 :
49  interfaceCompositionModel(dict, interface),
50  YNonVapour_
51  (
52  IOobject
53  (
54  IOobject::groupName("YNonVapour", this->interface().name()),
55  interface.mesh().time().name(),
56  interface.mesh()
57  ),
58  interface.mesh(),
60  ),
61  YNonVapourPrime_
62  (
63  IOobject
64  (
65  IOobject::groupName("YNonVapourPrime", this->interface().name()),
66  interface.mesh().time().name(),
67  interface.mesh()
68  ),
69  interface.mesh(),
71  )
72 {
74  {
75  speciesModels_.insert
76  (
77  *iter,
79  (
81  (
82  dict.subDict(*iter),
83  interface
84  )
85  )
86  );
87  }
88 }
89 
90 
91 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
92 
94 {}
95 
96 
97 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
98 
100 {
101  YNonVapour_ = scalar(1);
102 
103  forAllIter
104  (
106  speciesModels_,
107  iter
108  )
109  {
110  iter()->update(Tf);
111 
112  YNonVapour_ -=
113  otherComposition().Y(iter.key())
114  *iter()->Yf(iter.key(), Tf);
115 
116  YNonVapourPrime_ -=
117  otherComposition().Y(iter.key())
118  *iter()->YfPrime(iter.key(), Tf);
119  }
120 }
121 
122 
124 (
125  const word& speciesName,
126  const volScalarField& Tf
127 ) const
128 {
129  if (species().found(speciesName))
130  {
131  return
132  otherComposition().Y(speciesName)
133  *speciesModels_[speciesName]->Yf(speciesName, Tf);
134  }
135  else
136  {
137  return composition().Y(speciesName)*YNonVapour_;
138  }
139 }
140 
141 
144 (
145  const word& speciesName,
146  const volScalarField& Tf
147 ) const
148 {
149  if (species().found(speciesName))
150  {
151  return
152  otherComposition().Y(speciesName)
153  *speciesModels_[speciesName]->YfPrime(speciesName, Tf);
154  }
155  else
156  {
157  return composition().Y(speciesName)*YNonVapourPrime_;
158  }
159 }
160 
161 
162 // ************************************************************************* //
bool found
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
An STL-conforming hash table.
Definition: HashTable.H:127
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
A wordList with hashed indices for faster lookup by name.
Generic base class for interface composition models. These models describe the composition in phase 1...
const sidedPhaseInterface & interface() const
Return the interface.
static autoPtr< interfaceCompositionModel > New(const dictionary &dict, const phaseInterface &interface, const bool outer=true)
const hashedWordList & species() const
Return the transferring species names.
Raoult's law of ideal mixing. A separate composition model is given for each species....
Definition: Raoult.H:55
Raoult(const dictionary &dict, const phaseInterface &interface)
Construct from a dictionary and an interface.
Definition: Raoult.C:44
virtual ~Raoult()
Destructor.
Definition: Raoult.C:93
virtual void update(const volScalarField &Tf)
Update the composition.
Definition: Raoult.C:99
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const
The interface species fraction.
Definition: Raoult.C:124
virtual tmp< volScalarField > YfPrime(const word &speciesName, const volScalarField &Tf) const
The interface species fraction derivative w.r.t. temperature.
Definition: Raoult.C:144
Class to represent an interface between phases. Derivations can further specify the configuration of ...
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
addToRunTimeSelectionTable(interfaceCompositionModel, Henry, dictionary)
Namespace for OpenFOAM.
const dimensionSet dimless
const dimensionSet dimTemperature
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
dictionary dict
basicSpecieMixture & composition