semiPermeableBaffleMassFractionFvPatchScalarField.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) 2017-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 "fvPatchFieldMapper.H"
28 #include "volFields.H"
29 #include "surfaceFields.H"
31 #include "basicSpecieMixture.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
38 (
39  const fvPatch& p,
41 )
42 :
45 {}
46 
47 
50 (
51  const fvPatch& p,
53  const dictionary& dict
54 )
55 :
56  mappedPatchBase(p.patch(), NEARESTPATCHFACE, dict),
58 {}
59 
60 
63 (
65  const fvPatch& p,
67  const fvPatchFieldMapper& mapper
68 )
69 :
70  mappedPatchBase(p.patch(), ptf),
72 {}
73 
74 
77 (
80 )
81 :
82  mappedPatchBase(ptf.patch().patch(), ptf),
84 {}
85 
86 
87 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
88 
91 {
92  if (c_ == scalar(0))
93  {
94  return tmp<scalarField>(new scalarField(patch().size(), Zero));
95  }
96 
97  const word& YName = internalField().name();
98 
99  const fvPatch& nbrPatch = patch().boundaryMesh()[samplePolyPatch().index()];
100 
101  const fluidThermo& thermo =
102  db().lookupObject<fluidThermo>(fluidThermo::dictName);
103 
104  // Get the cell-mass fractions
105  const scalarField Yc(patchInternalField());
106  scalarField nbrYc
107  (
108  nbrPatch.lookupPatchField<volScalarField, scalar>(YName)
109  .patchInternalField()
110  );
112 
113  // Get the patch delta coefficients multiplied by the diffusivity
114  const thermophysicalTransportModel& ttm =
115  db().lookupObject<thermophysicalTransportModel>
116  (
117  thermophysicalTransportModel::typeName
118  );
119  const scalarField alphaEffDeltap
120  (
121  ttm.alphaEff(patch().index())*patch().deltaCoeffs()
122  );
123  scalarField nbrAlphaEffDeltap
124  (
125  ttm.alphaEff(nbrPatch.index())*nbrPatch.deltaCoeffs()
126  );
127  mappedPatchBase::map().distribute(nbrAlphaEffDeltap);
128 
129  // Get the specie molecular weight, if needed
130  scalar Wi = NaN;
131  if (property_ != massFraction)
132  {
133  const basicSpecieMixture& mixture = composition(db());
134  Wi = mixture.Wi(mixture.species()[YName]);
135  }
136 
137  // Get the mixture molecular weights, if needed
138  tmp<scalarField> tW, tNbrW;
140  {
141  tW = thermo.W(patch().index());
142  tNbrW = thermo.W(nbrPatch.index());
144  }
145 
146  // Construct coefficients that convert mass fraction to the property that
147  // drives the transfer
148  scalarField k(patch().size(), 1), nbrK(patch().size(), 1);
149  switch(property_)
150  {
151  case massFraction:
152  break;
153 
154  case moleFraction:
155  k *= tW/Wi;
156  nbrK *= tNbrW/Wi;
157  break;
158 
159  case molarConcentration:
160  {
161  k *= thermo.rho(patch().index())/Wi;
162  scalarField nbrRhop(thermo.rho(nbrPatch.index()));
163  mappedPatchBase::map().distribute(nbrRhop);
164  nbrK *= nbrRhop/Wi;
165  }
166  break;
167 
168  case partialPressure:
169  {
170  k *= thermo.p().boundaryField()[patch().index()]*tW/Wi;
171  scalarField nbrPp(thermo.p().boundaryField()[nbrPatch.index()]);
173  nbrK *= nbrPp*tNbrW/Wi;
174  }
175  break;
176  }
177 
178  // The transport is limited by both the difference in the amount of the
179  // specie across the baffle and the rates of diffusion on either side. The
180  // coefficients associated with these three transfer process are
181  // harmonically averaged to represent this limiting.
182  return
183  patch().magSf()
184  /(1/c_ + k/alphaEffDeltap + nbrK/nbrAlphaEffDeltap)
185  *(k*Yc - nbrK*nbrYc);
186 }
187 
188 
190 (
191  Ostream& os
192 ) const
193 {
196  writeEntry(os, "value", *this);
197 }
198 
199 
200 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
201 
202 namespace Foam
203 {
205  (
208  );
209 }
210 
211 // ************************************************************************* //
Foam::surfaceFields.
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:180
dictionary dict
fluidReactionThermo & thermo
Definition: createFields.H:28
virtual tmp< volScalarField > rho() const =0
Density [kg/m^3].
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:174
semiPermeableBaffleMassFractionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
virtual void write(Ostream &) const
Write as a dictionary.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:62
This is a mass-fraction boundary condition for a semi-permeable baffle.
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...
const mapDistribute & map() const
Return reference to the parallel distribution map.
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.
virtual tmp< scalarField > calcPhiYp() const
Return the flux of this species through the baffle.
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:138
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:53
Foam::fvPatchFieldMapper.
Determines a mapping between patch face centres and mesh cell or face centres and processors they&#39;re ...
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
static const word dictName
Name of the thermophysical properties dictionary.
Definition: basicThermo.H:129
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
virtual tmp< volScalarField > W() const =0
Molecular weight [kg/kmol].
const GeometricField::Patch & lookupPatchField(const word &name) const
Lookup and return the patchField of the named field from the.
Abstract base class for thermophysical transport models (RAS, LES and laminar).
Abstract base class for specie-transferring mass fraction boundary conditions.
const scalarField & deltaCoeffs() const
Return the face - cell distance coeffient.
Definition: fvPatch.C:170
virtual tmp< volScalarField > alphaEff() const =0
Effective thermal turbulent diffusivity of mixture [kg/m/s].
label index() const
Return the index of this patch in the boundaryMesh.
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
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
const polyPatch & samplePolyPatch() const
Get the patch on the region.
Namespace for OpenFOAM.
const speciesTable & species() const
Return the table of species.