semiPermeableBaffleMassFractionFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2017 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
31 #include "fvPatchFieldMapper.H"
32 #include "volFields.H"
33 #include "surfaceFields.H"
34 #include "turbulenceModel.H"
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
40 (
41  const fvPatch& p,
43 )
44 :
46  mixedFvPatchScalarField(p, iF),
47  c_(0),
48  phiName_("phi")
49 {
50  refValue() = Zero;
51  refGrad() = Zero;
52  valueFraction() = Zero;
53 }
54 
55 
58 (
59  const fvPatch& p,
61  const dictionary& dict
62 )
63 :
64  mappedPatchBase(p.patch(), NEARESTPATCHFACE, dict),
65  mixedFvPatchScalarField(p, iF),
66  c_(dict.getOrDefault<scalar>("c", 0)),
67  phiName_(dict.getOrDefault<word>("phi", "phi"))
68 {
69  this->readValueEntry(dict, IOobjectOption::MUST_READ);
70 
71  refValue() = Zero;
72  refGrad() = Zero;
73  valueFraction() = Zero;
74 }
75 
76 
79 (
81  const fvPatch& p,
83  const fvPatchFieldMapper& mapper
84 )
85 :
86  mappedPatchBase(p.patch(), ptf),
87  mixedFvPatchScalarField(ptf, p, iF, mapper),
88  c_(ptf.c_),
89  phiName_(ptf.phiName_)
90 {}
91 
92 
95 (
97 )
98 :
99  mappedPatchBase(ptf.patch().patch(), ptf),
100  mixedFvPatchScalarField(ptf),
101  c_(ptf.c_),
102  phiName_(ptf.phiName_)
103 {}
104 
105 
108 (
111 )
112 :
113  mappedPatchBase(ptf.patch().patch(), ptf),
114  mixedFvPatchScalarField(ptf, iF),
115  c_(ptf.c_),
116  phiName_(ptf.phiName_)
117 {}
118 
119 
120 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
121 
124 {
125  if (c_ == scalar(0))
126  {
127  return tmp<scalarField>::New(patch().size(), Zero);
128  }
129 
130  const word& YName = internalField().name();
131 
132  const label nbrPatchi = samplePolyPatch().index();
133  const fvPatch& nbrPatch = patch().boundaryMesh()[nbrPatchi];
134 
135  const auto& nbrYp = nbrPatch.lookupPatchField<volScalarField>(YName);
136  scalarField nbrYc(nbrYp.patchInternalField());
138 
139  return c_*patch().magSf()*(patchInternalField() - nbrYc);
140 }
141 
142 
144 {
145  if (updated())
146  {
147  return;
148  }
149 
150  const scalarField& phip =
151  patch().lookupPatchField<surfaceScalarField>(phiName_);
152 
153  const turbulenceModel& turbModel =
154  db().lookupObject<turbulenceModel>
155  (
157  );
158  const scalarField muEffp(turbModel.muEff(patch().index()));
159  const scalarField AMuEffp(patch().magSf()*muEffp);
160 
161  valueFraction() = phip/(phip - patch().deltaCoeffs()*AMuEffp);
162  refGrad() = - phiY()/AMuEffp;
163 
164  mixedFvPatchScalarField::updateCoeffs();
165 }
166 
167 
169 (
170  Ostream& os
171 ) const
172 {
175  os.writeEntryIfDifferent<scalar>("c", 0, c_);
176  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
178 }
179 
180 
181 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
182 
183 namespace Foam
184 {
186  (
188  semiPermeableBaffleMassFractionFvPatchScalarField
189  );
190 }
191 
192 // ************************************************************************* //
Foam::surfaceFields.
dictionary dict
semiPermeableBaffleMassFractionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute List data using default commsType, default flip/negate operator.
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
This is a mass-fraction boundary condition for a semi-permeable baffle.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
const mapDistribute & map() const
Return reference to the parallel distribution map.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
fvPatchField< scalar > fvPatchScalarField
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
Determines a mapping between patch face centres and mesh cell or face centres and processors they&#39;re ...
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:336
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
OBJstream os(runTime.globalPath()/outputName)
virtual void write(Ostream &os) const
Write as a dictionary.
tmp< scalarField > phiY() const
Return the flux of this species through the baffle.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127