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-2024 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_
45  (
46  dict.lookupOrDefault<scalar>("massFluxFraction", dimless, 1.0)
47  )
48 {
49  refValue() = 1.0;
50  refGrad() = 0.0;
51  valueFraction() = 0.0;
52 
53  if (dict.found("value"))
54  {
56  (
57  Field<scalar>("value", iF.dimensions(), dict, p.size())
58  );
59  }
60  else
61  {
63  }
64 }
65 
68 (
70  const fvPatch& p,
72  const fieldMapper& mapper
73 )
74 :
75  mixedFvPatchField<scalar>(ptf, p, iF, mapper),
76  phiName_(ptf.phiName_),
77  rhoName_(ptf.rhoName_),
78  massFluxFraction_(ptf.massFluxFraction_)
79 {}
80 
81 
84 (
87 )
88 :
89  mixedFvPatchField<scalar>(tppsf, iF),
90  phiName_(tppsf.phiName_),
91  rhoName_(tppsf.rhoName_),
92  massFluxFraction_(tppsf.massFluxFraction_)
93 {}
94 
95 
96 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
97 
99 (
100  const fvPatchScalarField& ptf,
101  const fieldMapper& mapper
102 )
103 {
104  mixedFvPatchField<scalar>::map(ptf, mapper);
105 }
106 
107 
109 (
110  const fvPatchScalarField& ptf
111 )
112 {
114 }
115 
116 
118 {
119  if (this->updated())
120  {
121  return;
122  }
123 
124  const label patchi = patch().index();
125 
127  db().lookupType<fluidThermophysicalTransportModel>
128  (
129  internalField().group()
130  );
131 
132  const fvsPatchField<scalar>& phip =
133  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
134 
135  const scalarField alphap
136  (
137  ttm.kappaEff(patchi)/ttm.thermo().Cpv().boundaryField()[patchi]
138  );
139 
140  refValue() = massFluxFraction_;
141  refGrad() = 0.0;
142 
143  valueFraction() =
144  1.0
145  /
146  (
147  1.0 +
148  alphap*patch().deltaCoeffs()*patch().magSf()/max(mag(phip), small)
149  );
150 
152 
153  if (debug)
154  {
155  scalar phi = gSum(-phip*(*this));
156 
157  Info<< patch().boundaryMesh().mesh().name() << ':'
158  << patch().name() << ':'
159  << this->internalField().name() << " :"
160  << " mass flux[Kg/s]:" << phi
161  << endl;
162  }
163 }
164 
165 
167 write(Ostream& os) const
168 {
170  writeEntry(os, "phi", phiName_);
171  writeEntry(os, "rho", rhoName_);
172  writeEntry(os, "massFluxFraction", massFluxFraction_);
173  writeEntry(os, "value", *this);
174 }
175 
176 
177 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
178 
179 namespace Foam
180 {
182  (
185  );
186 
187 }
188 
189 // ************************************************************************* //
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 dimensionSet & dimensions() const
Return dimensions.
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:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class for fluid thermophysical transport models RAS, LES and laminar.
virtual const fluidThermo & thermo() const =0
Access function to fluid thermophysical properties.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:249
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:202
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:83
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 scalarField & valueFraction()
virtual void map(const fvPatchField< Type > &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
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 fieldMapper &)
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:257
const dimensionSet dimless
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.