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 "fieldMapper.H"
28 #include "volFields.H"
29 #include "surfaceFields.H"
32 #include "mappedFvPatchBaseBase.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 fieldMapper& 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 mapper and the neighbouring mesh and patch
94  const mappedFvPatchBaseBase& mapper =
96  const fvPatch& nbrPatch = mapper.nbrFvPatch();
97 
98  const fluidThermo& thermo =
99  db().lookupObject<fluidThermo>(physicalProperties::typeName);
100 
101  // Get the cell-mass fractions
102  const scalarField Yc(patchInternalField());
103  const scalarField nbrYc
104  (
105  mapper.fromNeighbour
106  (
107  nbrPatch.lookupPatchField<volScalarField, scalar>(YName)
108  .patchInternalField()
109  )
110  );
111 
112  // Get the patch delta coefficients multiplied by the diffusivity
114  db().lookupType<fluidThermophysicalTransportModel>();
115 
116  const scalarField alphaEffDeltap
117  (
118  ttm.kappaEff(patch().index())*patch().deltaCoeffs()
119  /ttm.thermo().Cp().boundaryField()[patch().index()]
120  );
121 
122  const scalarField nbrAlphaEffDeltap
123  (
124  mapper.fromNeighbour
125  (
126  ttm.kappaEff(nbrPatch.index())*nbrPatch.deltaCoeffs()
127  /ttm.thermo().Cp().boundaryField()[nbrPatch.index()]
128  )
129  );
130 
131  // Get the specie molecular weight, if needed
132  scalar Wi = NaN;
133  if (property_ != massFraction)
134  {
135  const fluidMulticomponentThermo& thermo = this->thermo(db());
136  Wi = thermo.Wi(thermo.species()[YName]).value();
137  }
138 
139  // Get the mixture molecular weights, if needed
140  tmp<scalarField> tW, tNbrW;
141  if (property_ == moleFraction || property_ == partialPressure)
142  {
143  tW = thermo.W(patch().index());
144  tNbrW = mapper.fromNeighbour(thermo.W(nbrPatch.index()));
145  }
146 
147  // Construct coefficients that convert mass fraction to the property that
148  // drives the transfer
149  scalarField k(patch().size(), 1), nbrK(patch().size(), 1);
150  switch(property_)
151  {
152  case massFraction:
153  break;
154 
155  case moleFraction:
156  k *= tW/Wi;
157  nbrK *= tNbrW/Wi;
158  break;
159 
160  case molarConcentration:
161  {
162  k *= thermo.rho(patch().index())/Wi;
163  tmp<scalarField> nbrRhop =
164  mapper.fromNeighbour(thermo.rho(nbrPatch.index()));
165  nbrK *= nbrRhop/Wi;
166  }
167  break;
168 
169  case partialPressure:
170  {
171  k *= thermo.p().boundaryField()[patch().index()]*tW/Wi;
172  tmp<scalarField> nbrPp =
173  mapper.fromNeighbour
174  (
175  thermo.p().boundaryField()[nbrPatch.index()]
176  );
177  nbrK *= nbrPp*tNbrW/Wi;
178  }
179  break;
180  }
181 
182  // The transport is limited by both the difference in the amount of the
183  // specie across the baffle and the rates of diffusion on either side. The
184  // coefficients associated with these three transfer process are
185  // harmonically averaged to represent this limiting.
186  return
187  patch().magSf()
188  /(1/c_ + k/alphaEffDeltap + nbrK/nbrAlphaEffDeltap)
189  *(k*Yc - nbrK*nbrYc);
190 }
191 
192 
194 (
195  Ostream& os
196 ) const
197 {
199  writeEntry(os, "value", *this);
200 }
201 
202 
203 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
204 
205 namespace Foam
206 {
208  (
211  );
212 }
213 
214 // ************************************************************************* //
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
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
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:156
virtual const scalarField & deltaCoeffs() const
Return the face - cell distance coefficient.
Definition: fvPatch.C:176
const GeometricField::Patch & lookupPatchField(const word &name) const
Lookup and return the patchField of the named field from the.
Base class for fv patches that provide mapping between two fv patches.
const fvPatch & nbrFvPatch() const
Get the patch to map from.
static const mappedFvPatchBaseBase & getMap(const fvPatch &patch)
Cast the given fvPatch to a mappedFvPatchBaseBase. Handle errors.
static void validateMapForField(const PatchField &field, const FieldType &iF, const dictionary &context, const label froms=from::any)
Validate that the map exists and is appropriate for the given.
virtual const speciesTable & species() const =0
The table of species.
virtual dimensionedScalar Wi(const label speciei) const =0
Molecular weight [kg/kmol].
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
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.