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-2020 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 "psiReactionThermo.H"
31 #include "rhoReactionThermo.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
38 (
39  const fvPatch& p,
41 )
42 :
44 {}
45 
46 
49 (
50  const fvPatch& p,
52  const dictionary& dict
53 )
54 :
56 {}
57 
58 
61 (
63  const fvPatch& p,
65  const fvPatchFieldMapper& mapper
66 )
67 :
69 {}
70 
71 
74 (
76 )
77 :
79 {}
80 
81 
84 (
87 )
88 :
90 {}
91 
92 
93 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94 
97 {
98  if (c_ == scalar(0))
99  {
100  return tmp<scalarField>(new scalarField(patch().size(), Zero));
101  }
102 
103  const word& YName = internalField().name();
104 
105  const fluidThermo& thermo =
106  db().lookupObject<fluidThermo>(fluidThermo::dictName);
107 
108  // Get the cell-mass fraction
109  const scalarField Yc(patchInternalField());
110 
111  // Get the patch delta coefficients multiplied by the diffusivity
112  const thermophysicalTransportModel& ttm =
113  db().lookupObject<thermophysicalTransportModel>
114  (
115  thermophysicalTransportModel::typeName
116  );
117  const scalarField alphaEffDeltap
118  (
119  ttm.alphaEff(patch().index())*patch().deltaCoeffs()
120  );
121 
122  // Get the specie molecular weight, if needed
123  scalar Wi = NaN;
124  if (property_ != massFraction)
125  {
126  const basicSpecieMixture& mixture = composition(db());
127  Wi = mixture.Wi(mixture.species()[YName]);
128  }
129 
130  // Get the mixture molecular weights, if needed
131  tmp<scalarField> tW;
133  {
134  tW = thermo.W(patch().index());
135  }
136 
137  // Construct coefficients that convert mass fraction to the property that
138  // drives the transfer
139  scalarField k(patch().size(), 1);
140  switch(property_)
141  {
142  case massFraction:
143  break;
144 
145  case moleFraction:
146  k *= tW/Wi;
147  break;
148 
149  case molarConcentration:
150  k *= thermo.rho(patch().index())/Wi;
151  break;
152 
153  case partialPressure:
154  k *= thermo.p().boundaryField()[patch().index()]*tW/Wi;
155  break;
156  }
157 
158  // The transport is limited by both the difference in the amount of the
159  // specie available and the rate of diffusion to the wall. The coefficients
160  // associated with these two transfer process are harmonically averaged to
161  // represent this limiting.
162  return
163  patch().magSf()
164  /(1/c_ + k/alphaEffDeltap)
165  *k*Yc;
166 }
167 
168 
170 (
171  Ostream& os
172 ) const
173 {
175  writeEntry(os, "value", *this);
176 }
177 
178 
179 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
180 
181 namespace Foam
182 {
184  (
187  );
188 }
189 
190 // ************************************************************************* //
Foam::surfaceFields.
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:158
const Boundary & boundaryField() const
Return const-reference to the boundary field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
virtual volScalarField & p()
Pressure [Pa].
Definition: fluidThermo.C:79
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.
rhoReactionThermo & thermo
Definition: createFields.H:28
Specialization 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.
Macros for easy insertion into run-time selection tables.
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:123
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:48
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.