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  (
55  "template<class Thermo, class OtherThermo> "
56  "Foam::interfaceCompositionModels::Henry<Thermo, OtherThermo>:: "
57  "Henry "
58  "( "
59  "const dictionary& dict, "
60  "const phasePair& pair "
61  ")"
62  ) << "Differing number of species and solubilities"
63  << exit(FatalError);
64  }
65 }
66 
67 
68 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
69 
70 template<class Thermo, class OtherThermo>
72 {}
73 
74 
75 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
76 
77 template<class Thermo, class OtherThermo>
79 (
80  const volScalarField& Tf
81 )
82 {
83  YSolvent_ = scalar(1);
84 
85  forAllConstIter(hashedWordList, this->speciesNames_, iter)
86  {
87  YSolvent_ -= Yf(*iter, Tf);
88  }
89 }
90 
91 
92 template<class Thermo, class OtherThermo>
95 (
96  const word& speciesName,
97  const volScalarField& Tf
98 ) const
99 {
100  if (this->speciesNames_.contains(speciesName))
101  {
102  const label index = this->speciesNames_[speciesName];
103 
104  return
105  k_[index]
106  *this->otherThermo_.composition().Y(speciesName)
107  *this->otherThermo_.rhoThermo::rho()
108  /this->thermo_.rhoThermo::rho();
109  }
110  else
111  {
112  return
113  YSolvent_
114  *this->thermo_.composition().Y(speciesName);
115  }
116 }
117 
118 
119 template<class Thermo, class OtherThermo>
122 (
123  const word& speciesName,
124  const volScalarField& Tf
125 ) const
126 {
127  return tmp<volScalarField>
128  (
129  new volScalarField
130  (
131  IOobject
132  (
133  IOobject::groupName("YfPrime", this->pair_.name()),
134  this->pair_.phase1().mesh().time().timeName(),
135  this->pair_.phase1().mesh()
136  ),
137  this->pair_.phase1().mesh(),
139  )
140  );
141 }
142 
143 
144 // ************************************************************************* //
phaseModel & phase1
Definition: createFields.H:12
virtual tmp< volScalarField > YfPrime(const word &speciesName, const volScalarField &Tf) const
The interface species fraction derivative w.r.t. temperature.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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
virtual void update(const volScalarField &Tf)
Update the composition.
dynamicFvMesh & mesh
Henry(const dictionary &dict, const phasePair &pair)
Construct from components.
static word groupName(Name name, const word &group)
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
stressControl lookup("compactNormalStress") >> compactNormalStress
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
error FatalError
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:131
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
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:118
word timeName
Definition: getTimeIndex.H:3
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const
The interface species fraction.