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 "psiReactionThermo.H"
32 #include "rhoReactionThermo.H"
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
39 (
40  const fvPatch& p,
42 )
43 :
46 {}
47 
48 
51 (
52  const fvPatch& p,
54  const dictionary& dict
55 )
56 :
57  mappedPatchBase(p.patch(), NEARESTPATCHFACE, dict),
59 {}
60 
61 
64 (
66  const fvPatch& p,
68  const fvPatchFieldMapper& mapper
69 )
70 :
71  mappedPatchBase(p.patch(), ptf),
73 {}
74 
75 
78 (
80 )
81 :
82  mappedPatchBase(ptf.patch().patch(), ptf),
84 {}
85 
86 
89 (
92 )
93 :
94  mappedPatchBase(ptf.patch().patch(), ptf),
96 {}
97 
98 
99 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
100 
103 {
104  if (c_ == scalar(0))
105  {
106  return tmp<scalarField>(new scalarField(patch().size(), Zero));
107  }
108 
109  const word& YName = internalField().name();
110 
111  const fvPatch& nbrPatch = patch().boundaryMesh()[samplePolyPatch().index()];
112 
113  const fluidThermo& thermo =
114  db().lookupObject<fluidThermo>(fluidThermo::dictName);
115 
116  // Get the cell-mass fractions
117  const scalarField Yc(patchInternalField());
118  scalarField nbrYc
119  (
120  nbrPatch.lookupPatchField<volScalarField, scalar>(YName)
121  .patchInternalField()
122  );
124 
125  // Get the patch delta coefficients multiplied by the diffusivity
126  const thermophysicalTransportModel& ttm =
127  db().lookupObject<thermophysicalTransportModel>
128  (
129  thermophysicalTransportModel::typeName
130  );
131  const scalarField alphaEffDeltap
132  (
133  ttm.alphaEff(patch().index())*patch().deltaCoeffs()
134  );
135  scalarField nbrAlphaEffDeltap
136  (
137  ttm.alphaEff(nbrPatch.index())*nbrPatch.deltaCoeffs()
138  );
139  mappedPatchBase::map().distribute(nbrAlphaEffDeltap);
140 
141  // Get the specie molecular weight, if needed
142  scalar Wi = NaN;
143  if (property_ != massFraction)
144  {
145  const basicSpecieMixture& mixture = composition(db());
146  Wi = mixture.Wi(mixture.species()[YName]);
147  }
148 
149  // Get the mixture molecular weights, if needed
150  tmp<scalarField> tW, tNbrW;
152  {
153  tW = thermo.W(patch().index());
154  tNbrW = thermo.W(nbrPatch.index());
156  }
157 
158  // Construct coefficients that convert mass fraction to the property that
159  // drives the transfer
160  scalarField k(patch().size(), 1), nbrK(patch().size(), 1);
161  switch(property_)
162  {
163  case massFraction:
164  break;
165 
166  case moleFraction:
167  k *= tW/Wi;
168  nbrK *= tNbrW/Wi;
169  break;
170 
171  case molarConcentration:
172  {
173  k *= thermo.rho(patch().index())/Wi;
174  scalarField nbrRhop(thermo.rho(nbrPatch.index()));
175  mappedPatchBase::map().distribute(nbrRhop);
176  nbrK *= nbrRhop/Wi;
177  }
178  break;
179 
180  case partialPressure:
181  {
182  k *= thermo.p().boundaryField()[patch().index()]*tW/Wi;
183  scalarField nbrPp(thermo.p().boundaryField()[nbrPatch.index()]);
185  nbrK *= nbrPp*tNbrW/Wi;
186  }
187  break;
188  }
189 
190  // The transport is limited by both the difference in the amount of the
191  // specie across the baffle and the rates of diffusion on either side. The
192  // coefficients associated with these three transfer process are
193  // harmonically averaged to represent this limiting.
194  return
195  patch().magSf()
196  /(1/c_ + k/alphaEffDeltap + nbrK/nbrAlphaEffDeltap)
197  *(k*Yc - nbrK*nbrYc);
198 }
199 
200 
202 (
203  Ostream& os
204 ) const
205 {
208  writeEntry(os, "value", *this);
209 }
210 
211 
212 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
213 
214 namespace Foam
215 {
217  (
220  );
221 }
222 
223 // ************************************************************************* //
Foam::surfaceFields.
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:179
dictionary dict
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:173
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:158
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:174
virtual void write(Ostream &) const
Write as a dictionary.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
virtual volScalarField & p()
Pressure [Pa].
Definition: fluidThermo.C:79
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
static const basicSpecieMixture & composition(const objectRegistry &db)
Access the composition for the given database.
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...
const mapDistribute & map() const
Return reference to the parallel distribution map.
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
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:137
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:48
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
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
virtual tmp< volScalarField > W() const =0
Molecular weight [kg/kmol].
Abstract base class for thermophysical transport models (RAS, LES and laminar).
Abstract base class for specie-transferring mass fraction boundary conditions.
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const Type *=nullptr) const
Lookup and return the patchField of the named field from the.
const scalarField & deltaCoeffs() const
Return the face - cell distance coeffient.
Definition: fvPatch.C:164
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.