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 {
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 fvPatchScalarField& nbrYp =
136  nbrPatch.lookupPatchField<volScalarField, scalar>(YName);
137  scalarField nbrYc(nbrYp.patchInternalField());
139 
140  return c_*patch().magSf()*(patchInternalField() - nbrYc);
141 }
142 
143 
145 {
146  if (updated())
147  {
148  return;
149  }
150 
151  const scalarField& phip =
152  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
153 
154  const turbulenceModel& turbModel =
155  db().lookupObject<turbulenceModel>
156  (
158  );
159  const scalarField muEffp(turbModel.muEff(patch().index()));
160  const scalarField AMuEffp(patch().magSf()*muEffp);
161 
162  valueFraction() = phip/(phip - patch().deltaCoeffs()*AMuEffp);
163  refGrad() = - phiY()/AMuEffp;
164 
165  mixedFvPatchScalarField::updateCoeffs();
166 }
167 
168 
170 (
171  Ostream& os
172 ) const
173 {
176  os.writeEntryIfDifferent<scalar>("c", 0, c_);
177  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
178  writeEntry("value", os);
179 }
180 
181 
182 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
183 
184 namespace Foam
185 {
187  (
189  semiPermeableBaffleMassFractionFvPatchScalarField
190  );
191 }
192 
193 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:47
Foam::fvPatchScalarField
fvPatchField< scalar > fvPatchScalarField
Definition: fvPatchFieldsFwd.H:33
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Definition: fvPatchField.C:377
Foam::fvPatchField::operator=
virtual void operator=(const UList< Type > &)
Foam::Ostream::writeEntryIfDifferent
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Definition: Ostream.H:244
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:46
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::semiPermeableBaffleMassFractionFvPatchScalarField
This is a mass-fraction boundary condition for a semi-permeable baffle.
Definition: semiPermeableBaffleMassFractionFvPatchScalarField.H:123
Foam::semiPermeableBaffleMassFractionFvPatchScalarField::write
virtual void write(Ostream &) const
Definition: semiPermeableBaffleMassFractionFvPatchScalarField.C:163
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::semiPermeableBaffleMassFractionFvPatchScalarField::phiY
tmp< scalarField > phiY() const
Definition: semiPermeableBaffleMassFractionFvPatchScalarField.C:116
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:57
Foam::Zero
static constexpr const zero Zero
Definition: zero.H:131
Foam::mappedPatchBase::map
const mapDistribute & map() const
Definition: mappedPatchBaseI.H:192
Foam::mappedPatchBase::samplePolyPatch
const polyPatch & samplePolyPatch() const
Definition: mappedPatchBase.C:1655
Foam::turbulenceModel::propertiesName
static const word propertiesName
Definition: turbulenceModel.H:96
surfaceFields.H
Foam::surfaceFields.
fvPatchFieldMapper.H
Foam::mappedPatchBase
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
Definition: mappedPatchBase.H:108
Foam::Field
Generic templated field type.
Definition: Field.H:59
Foam::compressible::turbulenceModel
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Definition: turbulentFluidThermoModel.H:59
Foam::semiPermeableBaffleMassFractionFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Definition: semiPermeableBaffleMassFractionFvPatchScalarField.C:137
semiPermeableBaffleMassFractionFvPatchScalarField.H
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::fvPatchField::patchInternalField
virtual tmp< Field< Type > > patchInternalField() const
Definition: fvPatchField.C:226
Foam::semiPermeableBaffleMassFractionFvPatchScalarField::semiPermeableBaffleMassFractionFvPatchScalarField
semiPermeableBaffleMassFractionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Definition: semiPermeableBaffleMassFractionFvPatchScalarField.C:33
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:53
Foam::mapDistribute::distribute
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Definition: mapDistributeTemplates.C:147
Foam::fvPatch::lookupPatchField
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const Type *=nullptr) const
Definition: fvPatchFvMeshTemplates.C:27
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:59
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:119
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Definition: atmBoundaryLayer.C:26
Foam::patchIdentifier::index
label index() const noexcept
Definition: patchIdentifier.H:143
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:50
Foam::foamVersion::patch
const std::string patch
Foam::mappedPatchBase::write
virtual void write(Ostream &os) const
Definition: mappedPatchBase.C:1935
Foam::tmp::New
static tmp< T > New(Args &&... args)
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:41
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:52
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
turbulenceModel.H
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:50