Henry.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 | Copyright (C) 2015 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 "Henry.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class Thermo, class OtherThermo>
32 (
33  const dictionary& dict,
34  const phasePair& pair
35 )
36 :
37  InterfaceCompositionModel<Thermo, OtherThermo>(dict, pair),
38  k_(dict.lookup("k")),
39  YSolvent_
40  (
41  IOobject
42  (
43  IOobject::groupName("YSolvent", pair.name()),
44  pair.phase1().mesh().time().timeName(),
45  pair.phase1().mesh()
46  ),
47  pair.phase1().mesh(),
48  dimensionedScalar("one", dimless, 1)
49  )
50 {
51  if (k_.size() != this->speciesNames_.size())
52  {
54  << "Differing number of species and solubilities"
55  << exit(FatalError);
56  }
57 }
58 
59 
60 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
61 
62 template<class Thermo, class OtherThermo>
64 {}
65 
66 
67 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
68 
69 template<class Thermo, class OtherThermo>
71 (
72  const volScalarField& Tf
73 )
74 {
75  YSolvent_ = scalar(1);
76 
77  forAllConstIter(hashedWordList, this->speciesNames_, iter)
78  {
79  YSolvent_ -= Yf(*iter, Tf);
80  }
81 }
82 
83 
84 template<class Thermo, class OtherThermo>
87 (
88  const word& speciesName,
89  const volScalarField& Tf
90 ) const
91 {
92  if (this->speciesNames_.contains(speciesName))
93  {
94  const label index = this->speciesNames_[speciesName];
95 
96  return
97  k_[index]
98  *this->otherThermo_.composition().Y(speciesName)
99  *this->otherThermo_.rhoThermo::rho()
100  /this->thermo_.rhoThermo::rho();
101  }
102  else
103  {
104  return
105  YSolvent_
106  *this->thermo_.composition().Y(speciesName);
107  }
108 }
109 
110 
111 template<class Thermo, class OtherThermo>
114 (
115  const word& speciesName,
116  const volScalarField& Tf
117 ) const
118 {
119  return tmp<volScalarField>
120  (
121  new volScalarField
122  (
123  IOobject
124  (
125  IOobject::groupName("YfPrime", this->pair_.name()),
126  this->pair_.phase1().mesh().time().timeName(),
127  this->pair_.phase1().mesh()
128  ),
129  this->pair_.phase1().mesh(),
131  )
132  );
133 }
134 
135 
136 // ************************************************************************* //
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const
The interface species fraction.
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
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
IOobject(const word &name, const fileName &instance, const objectRegistry &registry, readOption r=NO_READ, writeOption w=NO_WRITE, bool registerObject=true)
Construct from name, instance, registry, io options.
Definition: IOobject.C:114
virtual tmp< volScalarField > YfPrime(const word &speciesName, const volScalarField &Tf) const
The interface species fraction derivative w.r.t. temperature.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
static word groupName(Name name, const word &group)
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
Henry(const dictionary &dict, const phasePair &pair)
Construct from components.
word timeName
Definition: getTimeIndex.H:3
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
phaseModel & phase1
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual void update(const volScalarField &Tf)
Update the composition.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:54