adsorptionMassFractionFvPatchScalarField.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-2021 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 
27 #include "volFields.H"
28 #include "surfaceFields.H"
30 #include "basicSpecieMixture.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
37 (
38  const fvPatch& p,
40 )
41 :
43 {}
44 
45 
48 (
49  const fvPatch& p,
51  const dictionary& dict
52 )
53 :
55 {}
56 
57 
60 (
62  const fvPatch& p,
64  const fvPatchFieldMapper& mapper
65 )
66 :
68 {}
69 
70 
73 (
76 )
77 :
79 {}
80 
81 
82 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
83 
86 {
87  if (c_ == scalar(0))
88  {
89  return tmp<scalarField>(new scalarField(patch().size(), Zero));
90  }
91 
92  const word& YName = internalField().name();
93 
94  const fluidThermo& thermo =
95  db().lookupObject<fluidThermo>(physicalProperties::typeName);
96 
97  // Get the cell-mass fraction
98  const scalarField Yc(patchInternalField());
99 
100  // Get the patch delta coefficients multiplied by the diffusivity
101  const thermophysicalTransportModel& ttm =
102  db().lookupObject<thermophysicalTransportModel>
103  (
104  thermophysicalTransportModel::typeName
105  );
106  const scalarField alphaEffDeltap
107  (
108  ttm.alphaEff(patch().index())*patch().deltaCoeffs()
109  );
110 
111  // Get the specie molecular weight, if needed
112  scalar Wi = NaN;
113  if (property_ != massFraction)
114  {
115  const basicSpecieMixture& mixture = composition(db());
116  Wi = mixture.Wi(mixture.species()[YName]);
117  }
118 
119  // Get the mixture molecular weights, if needed
120  tmp<scalarField> tW;
122  {
123  tW = thermo.W(patch().index());
124  }
125 
126  // Construct coefficients that convert mass fraction to the property that
127  // drives the transfer
128  scalarField k(patch().size(), 1);
129  switch(property_)
130  {
131  case massFraction:
132  break;
133 
134  case moleFraction:
135  k *= tW/Wi;
136  break;
137 
138  case molarConcentration:
139  k *= thermo.rho(patch().index())/Wi;
140  break;
141 
142  case partialPressure:
143  k *= thermo.p().boundaryField()[patch().index()]*tW/Wi;
144  break;
145  }
146 
147  // The transport is limited by both the difference in the amount of the
148  // specie available and the rate of diffusion to the wall. The coefficients
149  // associated with these two transfer process are harmonically averaged to
150  // represent this limiting.
151  return
152  patch().magSf()
153  /(1/c_ + k/alphaEffDeltap)
154  *k*Yc;
155 }
156 
157 
159 (
160  Ostream& os
161 ) const
162 {
164  writeEntry(os, "value", *this);
165 }
166 
167 
168 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
169 
170 namespace Foam
171 {
173  (
176  );
177 }
178 
179 // ************************************************************************* //
Foam::surfaceFields.
fluidReactionThermo & thermo
Definition: createFields.H:28
virtual tmp< volScalarField > rho() const =0
Density [kg/m^3].
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
const Boundary & boundaryField() const
Return const-reference to the boundary field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual scalar Wi(const label speciei) const =0
Molecular weight of the given specie [kg/kmol].
label k
Boltzmann constant.
Specialisation of basicMixture for a mixture consisting of a number for molecular species...
static const basicSpecieMixture & composition(const objectRegistry &db)
Access the composition for the given database.
virtual volScalarField & p()=0
Pressure [Pa].
Macros for easy insertion into run-time selection tables.
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:53
Foam::fvPatchFieldMapper.
static const zero Zero
Definition: zero.H:97
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
virtual tmp< volScalarField > W() const =0
Molecular weight [kg/kmol].
This is a mass-fraction boundary condition for an adsorbing wall.
Abstract base class for thermophysical transport models (RAS, LES and laminar).
Abstract base class for specie-transferring mass fraction boundary conditions.
virtual tmp< volScalarField > alphaEff() const =0
Effective thermal turbulent diffusivity of mixture [kg/m/s].
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A class for managing temporary objects.
Definition: PtrList.H:53
adsorptionMassFractionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
tmp< scalarField > calcPhiYp() const
Return the flux of this species.
Namespace for OpenFOAM.
const speciesTable & species() const
Return the table of species.