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-2018 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 
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
31 #include "turbulenceModel.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
37 (
38  const fvPatch& p,
40 )
41 :
43  mixedFvPatchScalarField(p, iF),
44  c_(0),
45  phiName_("phi")
46 {
47  refValue() = Zero;
48  refGrad() = Zero;
49  valueFraction() = Zero;
50 }
51 
52 
55 (
56  const fvPatch& p,
58  const dictionary& dict
59 )
60 :
61  mappedPatchBase(p.patch(), NEARESTPATCHFACE, dict),
62  mixedFvPatchScalarField(p, iF),
63  c_(dict.lookupOrDefault<scalar>("c", scalar(0))),
64  phiName_(dict.lookupOrDefault<word>("phi", "phi"))
65 {
66  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
67 
68  refValue() = Zero;
69  refGrad() = Zero;
70  valueFraction() = Zero;
71 }
72 
73 
76 (
78  const fvPatch& p,
80  const fvPatchFieldMapper& mapper
81 )
82 :
83  mappedPatchBase(p.patch(), ptf),
84  mixedFvPatchScalarField(ptf, p, iF, mapper),
85  c_(ptf.c_),
86  phiName_(ptf.phiName_)
87 {}
88 
89 
92 (
94 )
95 :
96  mappedPatchBase(ptf.patch().patch(), ptf),
97  mixedFvPatchScalarField(ptf),
98  c_(ptf.c_),
99  phiName_(ptf.phiName_)
100 {}
101 
102 
105 (
108 )
109 :
110  mappedPatchBase(ptf.patch().patch(), ptf),
111  mixedFvPatchScalarField(ptf, iF),
112  c_(ptf.c_),
113  phiName_(ptf.phiName_)
114 {}
115 
116 
117 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
118 
121 {
122  if (c_ == scalar(0))
123  {
124  return tmp<scalarField>(new scalarField(patch().size(), Zero));
125  }
126 
127  const word& YName = internalField().name();
128 
129  const label nbrPatchi = samplePolyPatch().index();
130  const fvPatch& nbrPatch = patch().boundaryMesh()[nbrPatchi];
131 
132  const fvPatchScalarField& nbrYp =
133  nbrPatch.lookupPatchField<volScalarField, scalar>(YName);
134  scalarField nbrYc(nbrYp.patchInternalField());
136 
137  return c_*patch().magSf()*(patchInternalField() - nbrYc);
138 }
139 
140 
142 {
143  if (updated())
144  {
145  return;
146  }
147 
148  const scalarField& phip =
149  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
150 
151  const turbulenceModel& turbModel =
152  db().lookupObject<turbulenceModel>
153  (
155  );
156  const scalarField muEffp(turbModel.muEff(patch().index()));
157  const scalarField AMuEffp(patch().magSf()*muEffp);
158 
159  valueFraction() = phip/(phip - patch().deltaCoeffs()*AMuEffp);
160  refGrad() = - phiY()/AMuEffp;
161 
162  mixedFvPatchScalarField::updateCoeffs();
163 }
164 
165 
167 (
168  Ostream& os
169 ) const
170 {
173  writeEntryIfDifferent<scalar>(os, "c", scalar(0), c_);
174  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
175  writeEntry("value", os);
176 }
177 
178 
179 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
180 
181 namespace Foam
182 {
184  (
187  );
188 }
189 
190 // ************************************************************************* //
Foam::surfaceFields.
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:185
dictionary dict
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:59
semiPermeableBaffleMassFractionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
virtual void write(Ostream &) const
Write as a dictionary.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
This is a mass-fraction boundary condition for a semi-permeable baffle.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:362
const mapDistribute & map() const
Return reference to the parallel distribution map.
Abstract base class for turbulence models (RAS, LES and laminar).
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:143
Foam::fvPatchFieldMapper.
Determines a mapping between patch face centres and mesh cell or face centres and processors they&#39;re ...
static const zero Zero
Definition: zero.H:97
virtual label size() const
Return size.
Definition: fvPatch.H:161
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
volScalarField scalarField(fieldObject, mesh)
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:224
tmp< scalarField > phiY() const
Return the flux of this species through the baffle.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const Type *=nullptr) const
Lookup and return the patchField of the named field from the.
label index() const
Return the index of this patch in the boundaryMesh.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
const polyPatch & samplePolyPatch() const
Get the patch on the region.
Namespace for OpenFOAM.