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-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 "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  const dictionary& dict
43 )
44 :
46 {
48  (
49  *this,
50  iF,
51  dict,
53  );
54 }
55 
56 
59 (
61  const fvPatch& p,
63  const fvPatchFieldMapper& mapper
64 )
65 :
67 {}
68 
69 
72 (
75 )
76 :
78 {}
79 
80 
81 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
82 
85 {
86  if (c_ == scalar(0))
87  {
88  return tmp<scalarField>(new scalarField(patch().size(), Zero));
89  }
90 
91  const word& YName = internalField().name();
92 
93  // Get the coupling information from the mappedPatchBase
94  const mappedPatchBase& mpp = mappedPatchBase::getMap(patch().patch());
95  const polyMesh& nbrMesh = mpp.nbrMesh();
96  const label nbrPatchi = mpp.nbrPolyPatch().index();
97  const fvPatch& nbrPatch =
98  refCast<const fvMesh>(nbrMesh).boundary()[nbrPatchi];
99 
100  const fluidThermo& thermo =
101  db().lookupObject<fluidThermo>(physicalProperties::typeName);
102 
103  // Get the cell-mass fractions
104  const scalarField Yc(patchInternalField());
105  const scalarField nbrYc
106  (
107  mpp.fromNeighbour
108  (
109  nbrPatch.lookupPatchField<volScalarField, scalar>(YName)
110  .patchInternalField()
111  )
112  );
113 
114  // Get the patch delta coefficients multiplied by the diffusivity
116  db().lookupType<fluidThermophysicalTransportModel>();
117 
118  const scalarField alphaEffDeltap
119  (
120  ttm.kappaEff(patch().index())*patch().deltaCoeffs()
121  /ttm.thermo().Cp().boundaryField()[patch().index()]
122  );
123 
124  const scalarField nbrAlphaEffDeltap
125  (
126  mpp.fromNeighbour
127  (
128  ttm.kappaEff(nbrPatch.index())*nbrPatch.deltaCoeffs()
129  /ttm.thermo().Cp().boundaryField()[nbrPatch.index()]
130  )
131  );
132 
133  // Get the specie molecular weight, if needed
134  scalar Wi = NaN;
135  if (property_ != massFraction)
136  {
137  const basicSpecieMixture& mixture = composition(db());
138  Wi = mixture.Wi(mixture.species()[YName]);
139  }
140 
141  // Get the mixture molecular weights, if needed
142  tmp<scalarField> tW, tNbrW;
143  if (property_ == moleFraction || property_ == partialPressure)
144  {
145  tW = thermo.W(patch().index());
146  tNbrW = mpp.fromNeighbour(thermo.W(nbrPatch.index()));
147  }
148 
149  // Construct coefficients that convert mass fraction to the property that
150  // drives the transfer
151  scalarField k(patch().size(), 1), nbrK(patch().size(), 1);
152  switch(property_)
153  {
154  case massFraction:
155  break;
156 
157  case moleFraction:
158  k *= tW/Wi;
159  nbrK *= tNbrW/Wi;
160  break;
161 
162  case molarConcentration:
163  {
164  k *= thermo.rho(patch().index())/Wi;
165  tmp<scalarField> nbrRhop =
166  mpp.fromNeighbour(thermo.rho(nbrPatch.index()));
167  nbrK *= nbrRhop/Wi;
168  }
169  break;
170 
171  case partialPressure:
172  {
173  k *= thermo.p().boundaryField()[patch().index()]*tW/Wi;
174  tmp<scalarField> nbrPp =
175  mpp.fromNeighbour
176  (
177  thermo.p().boundaryField()[nbrPatch.index()]
178  );
179  nbrK *= nbrPp*tNbrW/Wi;
180  }
181  break;
182  }
183 
184  // The transport is limited by both the difference in the amount of the
185  // specie across the baffle and the rates of diffusion on either side. The
186  // coefficients associated with these three transfer process are
187  // harmonically averaged to represent this limiting.
188  return
189  patch().magSf()
190  /(1/c_ + k/alphaEffDeltap + nbrK/nbrAlphaEffDeltap)
191  *(k*Yc - nbrK*nbrYc);
192 }
193 
194 
196 (
197  Ostream& os
198 ) const
199 {
201  writeEntry(os, "value", *this);
202 }
203 
204 
205 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
206 
207 namespace Foam
208 {
210  (
213  );
214 }
215 
216 // ************************************************************************* //
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...
Generic GeometricField class.
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
Specialisation of basicMixture for a mixture consisting of a number for molecular species.
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:160
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:57
virtual tmp< volScalarField > W() const =0
Molecular weight [kg/kmol].
virtual volScalarField & p()=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.
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:175
virtual const scalarField & deltaCoeffs() const
Return the face - cell distance coefficient.
Definition: fvPatch.C:170
const GeometricField::Patch & lookupPatchField(const word &name) const
Lookup and return the patchField of the named field from the.
Engine which provides mapping between two patches.
const polyPatch & nbrPolyPatch() const
Get the patch to map from.
static const mappedPatchBase & getMap(const polyPatch &patch)
Cast the given polyPatch to a mappedPatchBase. Handle errors.
const polyMesh & nbrMesh() const
Get the mesh for the region to map from.
static void validateMapForField(const PatchFieldType &field, const FieldType &iF, const dictionary &context, const label froms=from::any)
Validate that the map exists and is appropriate for the given.
tmp< Field< Type > > fromNeighbour(const Field< Type > &nbrFld) const
Map/interpolate the neighbour patch field to this patch.
label index() const
Return the index of this patch in the boundaryMesh.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
This is a mass-fraction boundary condition for a semi-permeable baffle.
semiPermeableBaffleMassFractionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual tmp< scalarField > calcPhiYp() const
Return the flux of this species through the baffle.
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
static Type NaN()
Return a primitive with all components set to NaN.
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
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
basicSpecieMixture & composition
volScalarField & p
fluidMulticomponentThermo & thermo
Definition: createFields.H:31
static const label differentPatch
Foam::surfaceFields.