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-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 "saturated.H"
27 #include "phasePair.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace interfaceCompositionModels
35 {
36  defineTypeNameAndDebug(saturated, 0);
38  (
39  interfaceCompositionModel,
40  saturated,
41  dictionary
42  );
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
51 {
52  const dimensionedScalar Wi
53  (
54  "W",
56  composition().Wi(saturatedIndex_)
57  );
58 
59  return Wi/thermo().W()/thermo().p();
60 }
61 
62 
63 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
64 
66 (
67  const dictionary& dict,
68  const phasePair& pair
69 )
70 :
71  interfaceCompositionModel(dict, pair),
72  saturatedName_(species()[0]),
73  saturatedIndex_(composition().species()[saturatedName_]),
74  saturationModel_
75  (
76  saturationModel::New
77  (
78  dict.subDict("saturationPressure"),
79  pair
80  )
81  )
82 {
83  if (species().size() != 1)
84  {
86  << "saturated model is suitable for one species only."
87  << exit(FatalError);
88  }
89 }
90 
91 
92 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
93 
95 {}
96 
97 
98 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
99 
101 (
102  const volScalarField& Tf
103 )
104 {}
105 
106 
108 (
109  const word& speciesName,
110  const volScalarField& Tf
111 ) const
112 {
113  if (saturatedName_ == speciesName)
114  {
115  return wRatioByP()*saturationModel_->pSat(Tf);
116  }
117  else
118  {
119  const label speciesIndex = composition().species()[speciesName];
120 
121  return
122  composition().Y()[speciesIndex]
123  *(scalar(1) - wRatioByP()*saturationModel_->pSat(Tf))
124  /max(scalar(1) - composition().Y()[saturatedIndex_], small);
125  }
126 }
127 
128 
131 (
132  const word& speciesName,
133  const volScalarField& Tf
134 ) const
135 {
136  if (saturatedName_ == speciesName)
137  {
138  return wRatioByP()*saturationModel_->pSatPrime(Tf);
139  }
140  else
141  {
142  const label speciesIndex = composition().species()[speciesName];
143 
144  return
145  - composition().Y()[speciesIndex]
146  *wRatioByP()*saturationModel_->pSatPrime(Tf)
147  /max(scalar(1) - composition().Y()[saturatedIndex_], small);
148  }
149 }
150 
151 
152 // ************************************************************************* //
const hashedWordList & species() const
Return the transferring species names.
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
basicSpecieMixture & composition
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const
The interface species fraction.
rhoReactionThermo & thermo
Definition: createFields.H:28
Macros for easy insertion into run-time selection tables.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
saturated(const dictionary &dict, const phasePair &pair)
Construct from components.
virtual void update(const volScalarField &Tf)
Update the composition.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
Definition: dimensionSets.H:53
defineTypeNameAndDebug(combustionModel, 0)
virtual tmp< volScalarField > YfPrime(const word &speciesName, const volScalarField &Tf) const
The interface species fraction derivative w.r.t. temperature.
const basicSpecieMixture & composition() const
Return the composition.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
PtrList< volScalarField > & Y
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: PtrList.H:53
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
tmp< volScalarField > wRatioByP() const
Constant of proportionality between partial pressure and mass.
Namespace for OpenFOAM.
const speciesTable & species() const
Return the table of species.