totalFlowRateAdvectiveDiffusiveFvPatchScalarField.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) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2018-2020,2022 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"
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
40 (
41  const fvPatch& p,
43 )
44 :
45  mixedFvPatchField<scalar>(p, iF),
46  phiName_("phi"),
47  rhoName_("none"),
48  massFluxFraction_(1.0)
49 {
50  refValue() = 0.0;
51  refGrad() = 0.0;
52  valueFraction() = 0.0;
53 }
54 
55 
58 (
59  const fvPatch& p,
61  const dictionary& dict
62 )
63 :
64  mixedFvPatchField<scalar>(p, iF),
65  phiName_(dict.getOrDefault<word>("phi", "phi")),
66  rhoName_(dict.getOrDefault<word>("rho", "none")),
67  massFluxFraction_(dict.getOrDefault<scalar>("massFluxFraction", 1))
68 {
69 
70  refValue() = 1.0;
71  refGrad() = 0.0;
72  valueFraction() = 0.0;
73 
74  if (!this->readValueEntry(dict))
75  {
77  }
78 }
79 
80 
83 (
84  const totalFlowRateAdvectiveDiffusiveFvPatchScalarField& ptf,
85  const fvPatch& p,
86  const DimensionedField<scalar, volMesh>& iF,
87  const fvPatchFieldMapper& mapper
88 )
89 :
90  mixedFvPatchField<scalar>(ptf, p, iF, mapper),
91  phiName_(ptf.phiName_),
92  rhoName_(ptf.rhoName_),
93  massFluxFraction_(ptf.massFluxFraction_)
94 {}
95 
96 
99 (
101 )
102 :
103  mixedFvPatchField<scalar>(tppsf),
104  phiName_(tppsf.phiName_),
105  rhoName_(tppsf.rhoName_),
106  massFluxFraction_(tppsf.massFluxFraction_)
107 {}
108 
109 
112 (
115 )
116 :
117  mixedFvPatchField<scalar>(tppsf, iF),
118  phiName_(tppsf.phiName_),
119  rhoName_(tppsf.rhoName_),
120  massFluxFraction_(tppsf.massFluxFraction_)
121 {}
122 
123 
124 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
125 
127 (
128  const fvPatchFieldMapper& m
129 )
130 {
132 }
133 
134 
136 (
137  const fvPatchScalarField& ptf,
138  const labelList& addr
139 )
140 {
142 }
143 
144 
146 {
147  if (this->updated())
148  {
149  return;
150  }
151 
152  const label patchi = patch().index();
153 
154  const auto& turbModel =
155  db().lookupObject<compressible::turbulenceModel>
156  (
158  (
160  internalField().group()
161  )
162  );
163 
164  const auto& phip = patch().lookupPatchField<surfaceScalarField>(phiName_);
165 
166  const scalarField alphap(turbModel.alphaEff(patchi));
167 
168  refValue() = massFluxFraction_;
169  refGrad() = 0.0;
170 
171  valueFraction() =
172  1.0
173  /(
174  1.0
175  + alphap*patch().deltaCoeffs()*patch().magSf()/max(mag(phip), SMALL)
176  );
177 
179 
180  if (debug)
181  {
182  scalar phi = gSum(-phip*(*this));
183 
184  Info<< patch().boundaryMesh().mesh().name() << ':'
185  << patch().name() << ':'
186  << this->internalField().name() << " :"
187  << " mass flux[Kg/s]:" << phi
188  << endl;
189  }
190 }
191 
192 
194 (
195  Ostream& os
196 ) const
197 {
199  os.writeEntry("phi", phiName_);
200  os.writeEntry("rho", rhoName_);
201  os.writeEntry("massFluxFraction", massFluxFraction_);
203 }
204 
205 
206 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
207 
208 namespace Foam
209 {
211  (
213  totalFlowRateAdvectiveDiffusiveFvPatchScalarField
214  );
215 
216 }
217 
218 // ************************************************************************* //
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the "value" entry into *this.
Definition: fvPatchField.C:32
Foam::surfaceFields.
This BC is used for species inlets. The diffusion and advection fluxes are considered to calculate th...
dictionary dict
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
This boundary condition provides a base class for &#39;mixed&#39; type boundary conditions, i.e. conditions that mix fixed value and patch-normal gradient conditions.
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
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
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Macros for easy insertion into run-time selection tables.
virtual Field< scalar > & refValue()
constexpr const char *const group
Group name for atomic constants.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Type gSum(const FieldField< Field, Type > &f)
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.
virtual scalarField & valueFraction()
A FieldMapper for finite-volume patch fields.
int debug
Static debugging option.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
OBJstream os(runTime.globalPath()/outputName)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:309
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:391
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.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
totalFlowRateAdvectiveDiffusiveFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
virtual Field< scalar > & refGrad()
Namespace for OpenFOAM.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.