specieTransferMassFractionFvPatchScalarField.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) 2019-2024 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 
28 #include "volFields.H"
29 #include "surfaceFields.H"
32 
33 // * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  template<>
38  const char* NamedEnum
39  <
41  4
42  >::names[] =
43  {
44  "massFraction",
45  "moleFraction",
46  "molarConcentration",
47  "partialPressure"
48  };
49 }
50 
51 const Foam::NamedEnum
52 <
54  4
56 
57 
58 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
59 
62 (
63  const objectRegistry& db
64 )
65 {
66  return
68  (
69  physicalProperties::typeName
70  );
71 }
72 
73 
74 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
75 
78 (
79  const fvPatch& p,
81  const dictionary& dict
82 )
83 :
84  mixedFvPatchScalarField(p, iF, dict, false),
85  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
86  UName_(dict.lookupOrDefault<word>("U", "U")),
87  phiYp_(p.size(), 0),
88  timeIndex_(-1),
89  c_(dict.lookupOrDefault<scalar>("c", unitAny, scalar(0))),
90  property_
91  (
92  c_ == scalar(0)
93  ? massFraction
94  : propertyNames_.read(dict.lookup("property"))
95  )
96 {
98  (
99  scalarField("value", iF.dimensions(), dict, p.size())
100  );
101 
102  refValue() = Zero;
103  refGrad() = Zero;
104  valueFraction() = Zero;
105 }
106 
107 
110 (
112  const fvPatch& p,
114  const fieldMapper& mapper
115 )
116 :
117  mixedFvPatchScalarField(ptf, p, iF, mapper),
118  phiName_(ptf.phiName_),
119  UName_(ptf.UName_),
120  phiYp_(p.size(), 0),
121  timeIndex_(-1),
122  c_(ptf.c_),
123  property_(ptf.property_)
124 {}
125 
126 
129 (
132 )
133 :
134  mixedFvPatchScalarField(ptf, iF),
135  phiName_(ptf.phiName_),
136  UName_(ptf.UName_),
137  phiYp_(ptf.size(), 0),
138  timeIndex_(-1),
139  c_(ptf.c_),
140  property_(ptf.property_)
141 {}
142 
143 
144 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
145 
147 (
148  const fvPatchScalarField& ptf,
149  const fieldMapper& mapper
150 )
151 {
152  mixedFvPatchScalarField::map(ptf, mapper);
153 
155  refCast<const specieTransferMassFractionFvPatchScalarField>(ptf);
156 
157  mapper(phiYp_, tiptf.phiYp_);
158 }
159 
160 
162 (
163  const fvPatchScalarField& ptf
164 )
165 {
166  mixedFvPatchScalarField::reset(ptf);
167 
169  refCast<const specieTransferMassFractionFvPatchScalarField>(ptf);
170 
171  phiYp_.reset(tiptf.phiYp_);
172 }
173 
174 
175 const Foam::scalarField&
177 {
178  if (timeIndex_ != this->db().time().timeIndex())
179  {
180  timeIndex_ = this->db().time().timeIndex();
181 
182  phiYp_ = calcPhiYp();
183  }
184 
185  return phiYp_;
186 }
187 
188 
190 {
191  if (updated())
192  {
193  return;
194  }
195 
196  // Get the fluxes
197  const scalarField& phip =
198  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
199  const fvPatchVectorField& Up =
200  patch().lookupPatchField<volVectorField, vector>(UName_);
201  tmp<scalarField> uPhip =
202  refCast<const specieTransferVelocityFvPatchVectorField>(Up).phip();
203 
205  db().lookupType<fluidThermophysicalTransportModel>();
206 
207  // Get the diffusivity
208  const scalarField AAlphaEffp
209  (
210  patch().magSf()
211  *ttm.kappaEff(patch().index())
212  /ttm.thermo().Cp().boundaryField()[patch().index()]
213  );
214 
215  // Set the gradient and value so that the transport and diffusion combined
216  // result in the desired specie flux
217  valueFraction() = phip/(phip - patch().deltaCoeffs()*AAlphaEffp);
218  refValue() = *this;
219  refGrad() = phip*(*this - phiYp()/uPhip)/AAlphaEffp;
220 
221  mixedFvPatchScalarField::updateCoeffs();
222 }
223 
224 
226 (
227  Ostream& os
228 ) const
229 {
231  writeEntryIfDifferent<scalar>(os, "c", scalar(0), c_);
232  writeEntry(os, "property", propertyNames_[property_]);
233  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
234  writeEntryIfDifferent<word>(os, "U", "U", UName_);
235 }
236 
237 
238 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:54
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
virtual const volScalarField & Cp() const =0
Heat capacity at constant pressure [J/kg/K].
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Base-class for multi-component fluid thermodynamic properties.
Abstract base class for fluid thermophysical transport models RAS, LES and laminar.
virtual const fluidThermo & thermo() const =0
Access function to fluid thermophysical properties.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
Registry of regIOobjects.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
Abstract base class for specie-transferring mass fraction boundary conditions.
virtual const scalarField & phiYp() const
Return the flux of this species.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
static const fluidMulticomponentThermo & thermo(const objectRegistry &db)
Get the thermo from the given database.
static const NamedEnum< property, 4 > propertyNames_
Property type names.
specieTransferMassFractionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual tmp< volScalarField > kappaEff() const =0
Effective thermal turbulent conductivity.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
bool read(const char *, int32_t &)
Definition: int32IO.C:85
VolField< vector > volVectorField
Definition: volFieldsFwd.H:65
const unitConversion unitAny
SurfaceField< scalar > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
label timeIndex
Definition: getTimeIndex.H:4
dictionary dict
volScalarField & p
Foam::surfaceFields.