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-2022 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"
32 #include "mappedPatchBase.H"
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
39 (
40  const fvPatch& p,
42 )
43 :
45 {}
46 
47 
50 (
51  const fvPatch& p,
53  const dictionary& dict
54 )
55 :
57 {
58  if (!isA<mappedPatchBase>(patch().patch()))
59  {
61  << "' not type '" << mappedPatchBase::typeName << "'"
62  << "\n for patch " << p.name()
63  << " of field " << internalField().name()
64  << " in file " << internalField().objectPath()
65  << exit(FatalError);
66  }
67 }
68 
69 
72 (
74  const fvPatch& p,
76  const fvPatchFieldMapper& mapper
77 )
78 :
80 {}
81 
82 
85 (
88 )
89 :
91 {}
92 
93 
94 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95 
98 {
99  if (c_ == scalar(0))
100  {
101  return tmp<scalarField>(new scalarField(patch().size(), Zero));
102  }
103 
104  const word& YName = internalField().name();
105 
106  // Get the coupling information from the mappedPatchBase
107  const mappedPatchBase& mpp =
108  refCast<const mappedPatchBase>(patch().patch());
109  const polyMesh& nbrMesh = mpp.sampleMesh();
110  const label samplePatchi = mpp.samplePolyPatch().index();
111  const fvPatch& nbrPatch =
112  refCast<const fvMesh>(nbrMesh).boundary()[samplePatchi];
113 
114  const fluidThermo& thermo =
115  db().lookupObject<fluidThermo>(physicalProperties::typeName);
116 
117  // Get the cell-mass fractions
118  const scalarField Yc(patchInternalField());
119  scalarField nbrYc
120  (
121  nbrPatch.lookupPatchField<volScalarField, scalar>(YName)
122  .patchInternalField()
123  );
124  mpp.distribute(nbrYc);
125 
126  // Get the patch delta coefficients multiplied by the diffusivity
127  const thermophysicalTransportModel& ttm =
128  db().lookupObject<thermophysicalTransportModel>
129  (
130  thermophysicalTransportModel::typeName
131  );
132  const scalarField alphaEffDeltap
133  (
134  ttm.alphaEff(patch().index())*patch().deltaCoeffs()
135  );
136  scalarField nbrAlphaEffDeltap
137  (
138  ttm.alphaEff(nbrPatch.index())*nbrPatch.deltaCoeffs()
139  );
140  mpp.distribute(nbrAlphaEffDeltap);
141 
142  // Get the specie molecular weight, if needed
143  scalar Wi = NaN;
144  if (property_ != massFraction)
145  {
146  const basicSpecieMixture& mixture = composition(db());
147  Wi = mixture.Wi(mixture.species()[YName]);
148  }
149 
150  // Get the mixture molecular weights, if needed
151  tmp<scalarField> tW, tNbrW;
153  {
154  tW = thermo.W(patch().index());
155  tNbrW = thermo.W(nbrPatch.index());
156  mpp.distribute(tNbrW.ref());
157  }
158 
159  // Construct coefficients that convert mass fraction to the property that
160  // drives the transfer
161  scalarField k(patch().size(), 1), nbrK(patch().size(), 1);
162  switch(property_)
163  {
164  case massFraction:
165  break;
166 
167  case moleFraction:
168  k *= tW/Wi;
169  nbrK *= tNbrW/Wi;
170  break;
171 
172  case molarConcentration:
173  {
174  k *= thermo.rho(patch().index())/Wi;
175  scalarField nbrRhop(thermo.rho(nbrPatch.index()));
176  mpp.distribute(nbrRhop);
177  nbrK *= nbrRhop/Wi;
178  }
179  break;
180 
181  case partialPressure:
182  {
183  k *= thermo.p().boundaryField()[patch().index()]*tW/Wi;
184  scalarField nbrPp(thermo.p().boundaryField()[nbrPatch.index()]);
185  mpp.distribute(nbrPp);
186  nbrK *= nbrPp*tNbrW/Wi;
187  }
188  break;
189  }
190 
191  // The transport is limited by both the difference in the amount of the
192  // specie across the baffle and the rates of diffusion on either side. The
193  // coefficients associated with these three transfer process are
194  // harmonically averaged to represent this limiting.
195  return
196  patch().magSf()
197  /(1/c_ + k/alphaEffDeltap + nbrK/nbrAlphaEffDeltap)
198  *(k*Yc - nbrK*nbrYc);
199 }
200 
201 
203 (
204  Ostream& os
205 ) const
206 {
208  writeEntry(os, "value", *this);
209 }
210 
211 
212 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
213 
214 namespace Foam
215 {
217  (
220  );
221 }
222 
223 // ************************************************************************* //
Foam::surfaceFields.
fluidReactionThermo & thermo
Definition: createFields.H:28
virtual tmp< volScalarField > rho() const =0
Density [kg/m^3].
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:175
semiPermeableBaffleMassFractionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
const polyMesh & sampleMesh() const
Get the region mesh.
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
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
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...
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.
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
faceListList boundary(nPatches)
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].
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).
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
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...
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
virtual const word & name() const
Return name.
Definition: fvPatch.H:145
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.