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-2023 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  const dictionary& dict
39 )
40 :
41  mixedFvPatchField<scalar>(p, iF, dict, false),
42  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
43  rhoName_(dict.lookupOrDefault<word>("rho", "none")),
44  massFluxFraction_(dict.lookupOrDefault<scalar>("massFluxFraction", 1.0))
45 {
46 
47  refValue() = 1.0;
48  refGrad() = 0.0;
49  valueFraction() = 0.0;
50 
51  if (dict.found("value"))
52  {
54  (
55  Field<scalar>("value", dict, p.size())
56  );
57  }
58  else
59  {
61  }
62 }
63 
66 (
68  const fvPatch& p,
70  const fvPatchFieldMapper& mapper
71 )
72 :
73  mixedFvPatchField<scalar>(ptf, p, iF, mapper),
74  phiName_(ptf.phiName_),
75  rhoName_(ptf.rhoName_),
76  massFluxFraction_(ptf.massFluxFraction_)
77 {}
78 
79 
82 (
85 )
86 :
87  mixedFvPatchField<scalar>(tppsf, iF),
88  phiName_(tppsf.phiName_),
89  rhoName_(tppsf.rhoName_),
90  massFluxFraction_(tppsf.massFluxFraction_)
91 {}
92 
93 
94 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95 
97 (
98  const fvPatchScalarField& ptf,
99  const fvPatchFieldMapper& mapper
100 )
101 {
102  mixedFvPatchField<scalar>::map(ptf, mapper);
103 }
104 
105 
107 (
108  const fvPatchScalarField& ptf
109 )
110 {
112 }
113 
114 
116 {
117  if (this->updated())
118  {
119  return;
120  }
121 
122  const label patchi = patch().index();
123 
125  db().lookupType<fluidThermophysicalTransportModel>
126  (
127  internalField().group()
128  );
129 
130  const fvsPatchField<scalar>& phip =
131  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
132 
133  const scalarField alphap
134  (
135  ttm.kappaEff(patchi)/ttm.thermo().Cpv().boundaryField()[patchi]
136  );
137 
138  refValue() = massFluxFraction_;
139  refGrad() = 0.0;
140 
141  valueFraction() =
142  1.0
143  /
144  (
145  1.0 +
146  alphap*patch().deltaCoeffs()*patch().magSf()/max(mag(phip), small)
147  );
148 
150 
151  if (debug)
152  {
153  scalar phi = gSum(-phip*(*this));
154 
155  Info<< patch().boundaryMesh().mesh().name() << ':'
156  << patch().name() << ':'
157  << this->internalField().name() << " :"
158  << " mass flux[Kg/s]:" << phi
159  << endl;
160  }
161 }
162 
163 
165 write(Ostream& os) const
166 {
168  writeEntry(os, "phi", phiName_);
169  writeEntry(os, "rho", rhoName_);
170  writeEntry(os, "massFluxFraction", massFluxFraction_);
171  writeEntry(os, "value", *this);
172 }
173 
174 
175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176 
177 namespace Foam
178 {
180  (
183  );
184 
185 }
186 
187 // ************************************************************************* //
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...
const Boundary & boundaryField() const
Return const-reference to the boundary field.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
virtual const volScalarField & Cpv() const =0
Heat capacity at constant pressure/volume [J/kg/K].
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Abstract base class for fluid thermophysical transport models RAS, LES and laminar.
virtual const fluidThermo & thermo() const =0
Access function to fluid thermophysical properties.
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:251
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:204
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:82
This boundary condition provides a base class for 'mixed' type boundary conditions,...
virtual Field< scalar > & refValue()
virtual void reset(const fvPatchField< Type > &)
Reset the fvPatchField to the given fvPatchField.
virtual void map(const fvPatchField< Type > &, const fvPatchFieldMapper &)
Map the given fvPatchField onto this fvPatchField.
virtual scalarField & valueFraction()
virtual Field< scalar > & refGrad()
virtual tmp< volScalarField > kappaEff() const =0
Effective thermal turbulent conductivity.
This BC is used for species inlets. The diffusion and advection fluxes are considered to calculate th...
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void map(const fvPatchScalarField &, const fvPatchFieldMapper &)
Map the given fvPatchField onto this fvPatchField.
totalFlowRateAdvectiveDiffusiveFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
A class for handling words, derived from string.
Definition: word.H:62
label patchi
Namespace for OpenFOAM.
Type gSum(const FieldField< Field, Type > &f)
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
SurfaceField< scalar > surfaceScalarField
messageStream Info
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p
Foam::surfaceFields.