totalFlowRateAdvectiveDiffusiveFvPatchScalarField.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) 2011-2020 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 
27 #include "surfaceFields.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
35 (
36  const fvPatch& p,
38 )
39 :
41  phiName_("phi"),
42  rhoName_("none"),
43  massFluxFraction_(1.0)
44 {
45  refValue() = 0.0;
46  refGrad() = 0.0;
47  valueFraction() = 0.0;
48 }
49 
50 
53 (
54  const fvPatch& p,
56  const dictionary& dict
57 )
58 :
60  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
61  rhoName_(dict.lookupOrDefault<word>("rho", "none")),
62  massFluxFraction_(dict.lookupOrDefault<scalar>("massFluxFraction", 1.0))
63 {
64 
65  refValue() = 1.0;
66  refGrad() = 0.0;
67  valueFraction() = 0.0;
68 
69  if (dict.found("value"))
70  {
72  (
73  Field<scalar>("value", dict, p.size())
74  );
75  }
76  else
77  {
79  }
80 }
81 
84 (
86  const fvPatch& p,
88  const fvPatchFieldMapper& mapper
89 )
90 :
91  mixedFvPatchField<scalar>(ptf, p, iF, mapper),
92  phiName_(ptf.phiName_),
93  rhoName_(ptf.rhoName_),
94  massFluxFraction_(ptf.massFluxFraction_)
95 {}
96 
97 
100 (
102 )
103 :
105  phiName_(tppsf.phiName_),
106  rhoName_(tppsf.rhoName_),
107  massFluxFraction_(tppsf.massFluxFraction_)
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 {
131  m(*this, *this);
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 thermophysicalTransportModel& ttm =
156  (
158  (
159  thermophysicalTransportModel::typeName,
160  internalField().group()
161  )
162  );
163 
164  const fvsPatchField<scalar>& phip =
165  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
166 
167  const scalarField alphap(ttm.alphaEff(patchi));
168 
169  refValue() = massFluxFraction_;
170  refGrad() = 0.0;
171 
172  valueFraction() =
173  1.0
174  /
175  (
176  1.0 +
177  alphap*patch().deltaCoeffs()*patch().magSf()/max(mag(phip), small)
178  );
179 
181 
182  if (debug)
183  {
184  scalar phi = gSum(-phip*(*this));
185 
186  Info<< patch().boundaryMesh().mesh().name() << ':'
187  << patch().name() << ':'
188  << this->internalField().name() << " :"
189  << " mass flux[Kg/s]:" << phi
190  << endl;
191  }
192 }
193 
194 
196 write(Ostream& os) const
197 {
199  writeEntry(os, "phi", phiName_);
200  writeEntry(os, "rho", rhoName_);
201  writeEntry(os, "massFluxFraction", massFluxFraction_);
202  writeEntry(os, "value", *this);
203 }
204 
205 
206 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
207 
208 namespace Foam
209 {
211  (
214  );
215 
216 }
217 
218 // ************************************************************************* //
Foam::surfaceFields.
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:179
const char *const group
Group name for atomic constants.
This BC is used for species inlets. The diffusion and advection fluxes are considered to calculate th...
dictionary dict
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:667
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
const scalarField & magSf() const
Return face area magnitudes.
Definition: fvPatch.C:136
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:173
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:356
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:280
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
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.
phi
Definition: pEqn.H:104
virtual Field< scalar > & refValue()
Type gSum(const FieldField< Field, Type > &f)
const fvMesh & mesh() const
Return the mesh reference.
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
virtual scalarField & valueFraction()
Foam::fvPatchFieldMapper.
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:325
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
Abstract base class for thermophysical transport models (RAS, LES and laminar).
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
const word & name() const
Return reference to name.
Definition: fvMesh.H:253
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.
const scalarField & deltaCoeffs() const
Return the face - cell distance coeffient.
Definition: fvPatch.C:164
virtual tmp< volScalarField > alphaEff() const =0
Effective thermal turbulent diffusivity of mixture [kg/m/s].
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:229
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:167
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual const word & name() const
Return name.
Definition: fvPatch.H:143
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)
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
virtual Field< scalar > & refGrad()
Namespace for OpenFOAM.
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:332