saturated.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 "saturated.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace interfaceCompositionModels
34 {
37  (
39  saturated,
41  );
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
50 {
51  return thermo().Wi(saturatedIndex_)/thermo().W()/thermo().p();
52 }
53 
54 
55 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
56 
58 (
59  const dictionary& dict,
60  const phaseInterface& interface
61 )
62 :
63  interfaceCompositionModel(dict, interface),
64  saturatedName_(species()[0]),
65  saturatedIndex_(thermo().species()[saturatedName_]),
66  saturationModel_(saturationPressureModel::New("pSat", dict))
67 {
68  if (species().size() != 1)
69  {
71  << "saturated model is suitable for one species only."
72  << exit(FatalError);
73  }
74 }
75 
76 
77 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
78 
80 {}
81 
82 
83 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
84 
86 (
87  const volScalarField& Tf
88 )
89 {}
90 
91 
93 (
94  const word& speciesName,
95  const volScalarField& Tf
96 ) const
97 {
98  if (saturatedName_ == speciesName)
99  {
100  return wRatioByP()*saturationModel_->pSat(Tf);
101  }
102  else
103  {
104  const label speciesIndex = thermo().species()[speciesName];
105 
106  return
107  thermo().Y()[speciesIndex]
108  *(scalar(1) - wRatioByP()*saturationModel_->pSat(Tf))
109  /max(scalar(1) - thermo().Y()[saturatedIndex_], small);
110  }
111 }
112 
113 
116 (
117  const word& speciesName,
118  const volScalarField& Tf
119 ) const
120 {
121  if (saturatedName_ == speciesName)
122  {
123  return wRatioByP()*saturationModel_->pSatPrime(Tf);
124  }
125  else
126  {
127  const label speciesIndex = thermo().species()[speciesName];
128 
129  return
130  - thermo().Y()[speciesIndex]
131  *wRatioByP()*saturationModel_->pSatPrime(Tf)
132  /max(scalar(1) - thermo().Y()[saturatedIndex_], small);
133  }
134 }
135 
136 
137 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
virtual tmp< volScalarField > W() const =0
Molecular weight [kg/kmol].
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
virtual const volScalarField & p() const =0
Pressure [Pa].
Generic base class for interface composition models. These models describe the composition in phase 1...
const rhoFluidMulticomponentThermo & thermo() const
Return the thermo.
const hashedWordList & species() const
Return the transferring species names.
Model which uses a saturation pressure model for a single species to calculate the interface composit...
Definition: saturated.H:55
virtual void update(const volScalarField &Tf)
Update the composition.
Definition: saturated.C:86
tmp< volScalarField > wRatioByP() const
Constant of proportionality between partial pressure and mass.
Definition: saturated.C:49
label saturatedIndex_
saturated species index
Definition: saturated.H:64
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const
The interface species fraction.
Definition: saturated.C:93
virtual tmp< volScalarField > YfPrime(const word &speciesName, const volScalarField &Tf) const
The interface species fraction derivative w.r.t. temperature.
Definition: saturated.C:116
saturated(const dictionary &dict, const phaseInterface &interface)
Construct from a dictionary and an interface.
Definition: saturated.C:58
virtual dimensionedScalar Wi(const label speciei) const =0
Molecular weight [kg/kmol].
Class to represent an interface between phases. Derivations can further specify the configuration of ...
Model to describe the dependence of saturation pressure on temperature, and vice versa.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
addToRunTimeSelectionTable(interfaceCompositionModel, Henry, dictionary)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
error FatalError
dictionary dict
PtrList< volScalarField > & Y
fluidMulticomponentThermo & thermo
Definition: createFields.H:31