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-------------------------------------------------------------------------------
11License
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
38Foam::semiPermeableBaffleMassFractionFvPatchScalarField::
39semiPermeableBaffleMassFractionFvPatchScalarField
40(
41 const fvPatch& p,
43)
44:
45 mappedPatchBase(p.patch()),
46 mixedFvPatchScalarField(p, iF),
47 c_(0),
48 phiName_("phi")
49{
50 refValue() = Zero;
51 refGrad() = Zero;
52 valueFraction() = Zero;
53}
54
55
56Foam::semiPermeableBaffleMassFractionFvPatchScalarField::
57semiPermeableBaffleMassFractionFvPatchScalarField
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
77Foam::semiPermeableBaffleMassFractionFvPatchScalarField::
78semiPermeableBaffleMassFractionFvPatchScalarField
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
93Foam::semiPermeableBaffleMassFractionFvPatchScalarField::
94semiPermeableBaffleMassFractionFvPatchScalarField
95(
97)
98:
99 mappedPatchBase(ptf.patch().patch(), ptf),
100 mixedFvPatchScalarField(ptf),
101 c_(ptf.c_),
102 phiName_(ptf.phiName_)
103{}
104
105
106Foam::semiPermeableBaffleMassFractionFvPatchScalarField::
107semiPermeableBaffleMassFractionFvPatchScalarField
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
184namespace Foam
185{
187 (
190 );
191}
192
193// ************************************************************************* //
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...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
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:251
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:237
virtual void operator=(const UList< scalar > &)
Definition: fvPatchField.C:408
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const Type *=nullptr) const
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
const mapDistribute & map() const
Return reference to the parallel distribution map.
This is a mass-fraction boundary condition for a semi-permeable baffle.
tmp< scalarField > phiY() const
Return the flux of this species through the baffle.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A class for managing temporary objects.
Definition: tmp.H:65
Abstract base class for turbulence models (RAS, LES and laminar).
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dictionary dict
Foam::surfaceFields.