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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
30 #include "fvPatchFieldMapper.H"
31 #include "volFields.H"
32 #include "surfaceFields.H"
33 #include "turbulenceModel.H"
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
39 (
40  const fvPatch& p,
42 )
43 :
44  mappedPatchBase(p.patch()),
45  mixedFvPatchScalarField(p, iF),
46  c_(0),
47  phiName_("phi")
48 {
49  refValue() = Zero;
50  refGrad() = Zero;
51  valueFraction() = Zero;
52 }
53 
54 
57 (
58  const fvPatch& p,
60  const dictionary& dict
61 )
62 :
63  mappedPatchBase(p.patch(), NEARESTPATCHFACE, dict),
64  mixedFvPatchScalarField(p, iF),
65  c_(dict.lookupOrDefault<scalar>("c", 0)),
66  phiName_(dict.lookupOrDefault<word>("phi", "phi"))
67 {
68  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
69 
70  refValue() = Zero;
71  refGrad() = Zero;
72  valueFraction() = Zero;
73 }
74 
75 
78 (
80  const fvPatch& p,
82  const fvPatchFieldMapper& mapper
83 )
84 :
85  mappedPatchBase(p.patch(), ptf),
86  mixedFvPatchScalarField(ptf, p, iF, mapper),
87  c_(ptf.c_),
88  phiName_(ptf.phiName_)
89 {}
90 
91 
94 (
96 )
97 :
98  mappedPatchBase(ptf.patch().patch(), ptf),
99  mixedFvPatchScalarField(ptf),
100  c_(ptf.c_),
101  phiName_(ptf.phiName_)
102 {}
103 
104 
107 (
110 )
111 :
112  mappedPatchBase(ptf.patch().patch(), ptf),
113  mixedFvPatchScalarField(ptf, iF),
114  c_(ptf.c_),
115  phiName_(ptf.phiName_)
116 {}
117 
118 
119 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
120 
123 {
124  if (c_ == scalar(0))
125  {
126  return tmp<scalarField>::New(patch().size(), Zero);
127  }
128 
129  const word& YName = internalField().name();
130 
131  const label nbrPatchi = samplePolyPatch().index();
132  const fvPatch& nbrPatch = patch().boundaryMesh()[nbrPatchi];
133 
134  const fvPatchScalarField& nbrYp =
135  nbrPatch.lookupPatchField<volScalarField, scalar>(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, scalar>(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_);
177  writeEntry("value", os);
178 }
179 
180 
181 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
182 
183 namespace Foam
184 {
186  (
189  );
190 }
191 
192 // ************************************************************************* //
Foam::fvPatchField< scalar >
volFields.H
Foam::fvPatchField< scalar >::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:364
Foam::Ostream::writeEntryIfDifferent
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:231
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::mappedPatchBase::write
virtual void write(Ostream &) const
Write as a dictionary.
Definition: mappedPatchBase.C:1364
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::semiPermeableBaffleMassFractionFvPatchScalarField
This is a mass-fraction boundary condition for a semi-permeable baffle.
Definition: semiPermeableBaffleMassFractionFvPatchScalarField.H:127
Foam::semiPermeableBaffleMassFractionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: semiPermeableBaffleMassFractionFvPatchScalarField.C:169
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::semiPermeableBaffleMassFractionFvPatchScalarField::phiY
tmp< scalarField > phiY() const
Return the flux of this species through the baffle.
Definition: semiPermeableBaffleMassFractionFvPatchScalarField.C:122
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::mappedPatchBase::map
const mapDistribute & map() const
Return reference to the parallel distribution map.
Definition: mappedPatchBaseI.H:146
Foam::mappedPatchBase::samplePolyPatch
const polyPatch & samplePolyPatch() const
Get the patch on the region.
Definition: mappedPatchBase.C:1223
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
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:105
Foam::label
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:62
Foam::Field< scalar >
Foam::semiPermeableBaffleMassFractionFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: semiPermeableBaffleMassFractionFvPatchScalarField.C:143
semiPermeableBaffleMassFractionFvPatchScalarField.H
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Foam::fvPatchField::patchInternalField
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:225
Foam::semiPermeableBaffleMassFractionFvPatchScalarField::semiPermeableBaffleMassFractionFvPatchScalarField
semiPermeableBaffleMassFractionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: semiPermeableBaffleMassFractionFvPatchScalarField.C:39
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::mapDistribute::distribute
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
Definition: mapDistributeTemplates.C:152
Foam::fvPatch::lookupPatchField
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const Type *=nullptr) const
Definition: fvPatchFvMeshTemplates.C:34
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:63
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::patchIdentifier::index
label index() const
The index of this patch in the boundaryMesh.
Definition: patchIdentifier.H:133
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:54