Henry.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-2022 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"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace interfaceCompositionModels
34 {
35  defineTypeNameAndDebug(Henry, 0);
36  addToRunTimeSelectionTable(interfaceCompositionModel, Henry, dictionary);
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const dictionary& dict,
46  const phaseInterface& interface
47 )
48 :
49  interfaceCompositionModel(dict, interface),
50  k_(dict.lookup("k")),
51  YSolvent_
52  (
53  IOobject
54  (
55  IOobject::groupName("YSolvent", this->interface().name()),
56  interface.mesh().time().timeName(),
57  interface.mesh()
58  ),
59  interface.mesh(),
61  )
62 {
63  if (k_.size() != species().size())
64  {
66  << "Differing number of species and solubilities"
67  << exit(FatalError);
68  }
69 }
70 
71 
72 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
73 
75 {}
76 
77 
78 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
79 
81 {
82  YSolvent_ = scalar(1);
83 
84  forAllConstIter(hashedWordList, species(), iter)
85  {
86  YSolvent_ -= Yf(*iter, Tf);
87  }
88 }
89 
90 
92 (
93  const word& speciesName,
94  const volScalarField& Tf
95 ) const
96 {
97  if (species().found(speciesName))
98  {
99  const label index = species()[speciesName];
100 
101  return
102  k_[index]
103  *otherComposition().Y(speciesName)
104  *otherThermo().rho()
105  /thermo().rho();
106  }
107  else
108  {
109  return YSolvent_*composition().Y(speciesName);
110  }
111 }
112 
113 
115 (
116  const word& speciesName,
117  const volScalarField& Tf
118 ) const
119 {
120  return volScalarField::New
121  (
122  IOobject::groupName("YfPrime", interface().name()),
123  interface().mesh(),
125  );
126 }
127 
128 
129 // ************************************************************************* //
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
fluidReactionThermo & thermo
Definition: createFields.H:28
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:306
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
basicSpecieMixture & composition
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
Henry(const dictionary &dict, const phaseInterface &interface)
Construct from a dictionary and an interface.
const dimensionSet dimless
fvMesh & mesh
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const
The interface species fraction.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
stressControl lookup("compactNormalStress") >> compactNormalStress
static word groupName(Name name, const word &group)
word timeName
Definition: getTimeIndex.H:3
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual void update(const volScalarField &Tf)
Update the composition.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual tmp< volScalarField > YfPrime(const word &speciesName, const volScalarField &Tf) const
The interface species fraction derivative w.r.t. temperature.
A class for managing temporary objects.
Definition: PtrList.H:53
bool found
const dimensionSet dimTemperature
Namespace for OpenFOAM.