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-2025 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 const Foam::NamedEnum
36 <
38  4
40 {
41  "massFraction",
42  "moleFraction",
43  "molarConcentration",
44  "partialPressure"
45 };
46 
47 
48 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
49 
52 (
53  const objectRegistry& db
54 )
55 {
56  return
58  (
59  physicalProperties::typeName
60  );
61 }
62 
63 
64 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
65 
68 (
69  const fvPatch& p,
71  const dictionary& dict
72 )
73 :
74  mixedFvPatchScalarField(p, iF, dict, false),
75  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
76  UName_(dict.lookupOrDefault<word>("U", "U")),
77  phiYp_(p.size(), 0),
78  timeIndex_(-1),
79  c_(dict.lookupOrDefault<scalar>("c", unitAny, scalar(0))),
80  property_
81  (
82  c_ == scalar(0)
83  ? massFraction
84  : propertyNames_.read(dict.lookup("property"))
85  )
86 {
88  (
89  scalarField("value", iF.dimensions(), dict, p.size())
90  );
91 
92  refValue() = Zero;
93  refGrad() = Zero;
94  valueFraction() = Zero;
95 }
96 
97 
100 (
102  const fvPatch& p,
104  const fieldMapper& mapper
105 )
106 :
107  mixedFvPatchScalarField(ptf, p, iF, mapper),
108  phiName_(ptf.phiName_),
109  UName_(ptf.UName_),
110  phiYp_(p.size(), 0),
111  timeIndex_(-1),
112  c_(ptf.c_),
113  property_(ptf.property_)
114 {}
115 
116 
119 (
122 )
123 :
124  mixedFvPatchScalarField(ptf, iF),
125  phiName_(ptf.phiName_),
126  UName_(ptf.UName_),
127  phiYp_(ptf.size(), 0),
128  timeIndex_(-1),
129  c_(ptf.c_),
130  property_(ptf.property_)
131 {}
132 
133 
134 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
135 
137 (
138  const fvPatchScalarField& ptf,
139  const fieldMapper& mapper
140 )
141 {
142  mixedFvPatchScalarField::map(ptf, mapper);
143 
145  refCast<const specieTransferMassFractionFvPatchScalarField>(ptf);
146 
147  mapper(phiYp_, tiptf.phiYp_);
148 }
149 
150 
152 (
153  const fvPatchScalarField& ptf
154 )
155 {
156  mixedFvPatchScalarField::reset(ptf);
157 
159  refCast<const specieTransferMassFractionFvPatchScalarField>(ptf);
160 
161  phiYp_.reset(tiptf.phiYp_);
162 }
163 
164 
165 const Foam::scalarField&
167 {
168  if (timeIndex_ != this->db().time().timeIndex())
169  {
170  timeIndex_ = this->db().time().timeIndex();
171 
172  phiYp_ = calcPhiYp();
173  }
174 
175  return phiYp_;
176 }
177 
178 
180 {
181  if (updated())
182  {
183  return;
184  }
185 
186  // Get the fluxes
187  const scalarField& phip =
188  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
189  const fvPatchVectorField& Up =
190  patch().lookupPatchField<volVectorField, vector>(UName_);
191  tmp<scalarField> uPhip =
192  refCast<const specieTransferVelocityFvPatchVectorField>(Up).phip();
193 
195  db().lookupType<fluidThermophysicalTransportModel>();
196 
197  const volScalarField& Yi = refCast<const volScalarField>(internalField());
198 
199  // Get the diffusivity
200  const scalarField ADEffp
201  (
202  patch().magSf()*ttm.DEff(Yi, patch().index())
203  );
204 
205  // Compute the flux that we need to recover
206  const scalarField phiYp
207  (
208  this->phiYp()
209  - ADEffp*snGrad()
210  - patch().magSf()*ttm.j(Yi, patch().index())
211  );
212 
213  // Set the gradient and value so that the transport and diffusion combined
214  // result in the desired specie flux
215  valueFraction() = phip/(phip - patch().deltaCoeffs()*ADEffp);
216  refValue() = *this;
217  refGrad() = phip*(*this - phiYp/uPhip)/ADEffp;
218 
219  mixedFvPatchScalarField::updateCoeffs();
220 }
221 
222 
224 (
225  Ostream& os
226 ) const
227 {
229  writeEntryIfDifferent<scalar>(os, "c", scalar(0), c_);
230  writeEntry(os, "property", propertyNames_[property_]);
231  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
232  writeEntryIfDifferent<word>(os, "U", "U", UName_);
233 }
234 
235 
236 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:55
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
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 tmp< volScalarField > DEff(const volScalarField &Yi) const =0
Effective mass diffusion coefficient.
virtual tmp< surfaceScalarField > j(const volScalarField &Yi) const =0
Return the specie flux for the given specie mass-fraction [kg/m^2/s].
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:91
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:237
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.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
Definition: fvcSnGrad.C:45
static const zero Zero
Definition: zero.H:97
bool read(const char *, int32_t &)
Definition: int32IO.C:85
VolField< vector > volVectorField
Definition: volFieldsFwd.H:63
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.