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-2023 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"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
37 (
38  const fvPatch& p,
40  const dictionary& dict
41 )
42 :
44 {}
45 
46 
49 (
51  const fvPatch& p,
53  const fieldMapper& mapper
54 )
55 :
57 {}
58 
59 
62 (
65 )
66 :
68 {}
69 
70 
71 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
72 
75 {
76  if (c_ == scalar(0))
77  {
78  return tmp<scalarField>(new scalarField(patch().size(), Zero));
79  }
80 
81  const word& YName = internalField().name();
82 
83  const fluidThermo& thermo =
84  db().lookupObject<fluidThermo>(physicalProperties::typeName);
85 
86  // Get the cell-mass fraction
87  const scalarField Yc(patchInternalField());
88 
89  // Get the patch delta coefficients multiplied by the diffusivity
91  db().lookupType<fluidThermophysicalTransportModel>();
92 
93  const scalarField alphaEffDeltap
94  (
95  ttm.kappaEff(patch().index())*patch().deltaCoeffs()
96  /ttm.thermo().Cp().boundaryField()[patch().index()]
97  );
98 
99  // Get the specie molecular weight, if needed
100  scalar Wi = NaN;
101  if (property_ != massFraction)
102  {
103  const fluidMulticomponentThermo& thermo = this->thermo(db());
104  Wi = thermo.Wi(thermo.species()[YName]).value();
105  }
106 
107  // Get the mixture molecular weights, if needed
108  tmp<scalarField> tW;
109  if (property_ == moleFraction || property_ == partialPressure)
110  {
111  tW = thermo.W(patch().index());
112  }
113 
114  // Construct coefficients that convert mass fraction to the property that
115  // drives the transfer
116  scalarField k(patch().size(), 1);
117  switch(property_)
118  {
119  case massFraction:
120  break;
121 
122  case moleFraction:
123  k *= tW/Wi;
124  break;
125 
126  case molarConcentration:
127  k *= thermo.rho(patch().index())/Wi;
128  break;
129 
130  case partialPressure:
131  k *= thermo.p().boundaryField()[patch().index()]*tW/Wi;
132  break;
133  }
134 
135  // The transport is limited by both the difference in the amount of the
136  // specie available and the rate of diffusion to the wall. The coefficients
137  // associated with these two transfer process are harmonically averaged to
138  // represent this limiting.
139  return
140  patch().magSf()
141  /(1/c_ + k/alphaEffDeltap)
142  *k*Yc;
143 }
144 
145 
147 (
148  Ostream& os
149 ) const
150 {
152  writeEntry(os, "value", *this);
153 }
154 
155 
156 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
157 
158 namespace Foam
159 {
161  (
164  );
165 }
166 
167 // ************************************************************************* //
label k
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Boundary & boundaryField() const
Return const-reference to the boundary field.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
This is a mass-fraction boundary condition for an adsorbing wall.
tmp< scalarField > calcPhiYp() const
Return the flux of this species.
adsorptionMassFractionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual tmp< volScalarField > W() const =0
Molecular weight [kg/kmol].
virtual const volScalarField & Cp() const =0
Heat capacity at constant pressure [J/kg/K].
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:162
const Type & value() const
Return const reference to value.
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Base-class for multi-component fluid thermodynamic properties.
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:57
virtual const volScalarField & p() const =0
Pressure [Pa].
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
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
virtual const speciesTable & species() const =0
The table of species.
virtual dimensionedScalar Wi(const label speciei) const =0
Molecular weight [kg/kmol].
Abstract base class for specie-transferring mass fraction boundary conditions.
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
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p
fluidMulticomponentThermo & thermo
Definition: createFields.H:31
Foam::surfaceFields.