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-2022 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 (
103 )
104 :
105  mixedFvPatchField<scalar>(tppsf, iF),
106  phiName_(tppsf.phiName_),
107  rhoName_(tppsf.rhoName_),
108  massFluxFraction_(tppsf.massFluxFraction_)
109 {}
110 
111 
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113 
115 (
116  const fvPatchFieldMapper& m
117 )
118 {
119  m(*this, *this);
120 }
121 
122 
124 (
125  const fvPatchScalarField& ptf,
126  const labelList& addr
127 )
128 {
130 }
131 
132 
134 (
135  const fvPatchScalarField& ptf
136 )
137 {
139 }
140 
141 
143 {
144  if (this->updated())
145  {
146  return;
147  }
148 
149  const label patchi = patch().index();
150 
151  const thermophysicalTransportModel& ttm =
153  (
155  (
156  thermophysicalTransportModel::typeName,
157  internalField().group()
158  )
159  );
160 
161  const fvsPatchField<scalar>& phip =
162  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
163 
164  const scalarField alphap(ttm.alphaEff(patchi));
165 
166  refValue() = massFluxFraction_;
167  refGrad() = 0.0;
168 
169  valueFraction() =
170  1.0
171  /
172  (
173  1.0 +
174  alphap*patch().deltaCoeffs()*patch().magSf()/max(mag(phip), small)
175  );
176 
178 
179  if (debug)
180  {
181  scalar phi = gSum(-phip*(*this));
182 
183  Info<< patch().boundaryMesh().mesh().name() << ':'
184  << patch().name() << ':'
185  << this->internalField().name() << " :"
186  << " mass flux[Kg/s]:" << phi
187  << endl;
188  }
189 }
190 
191 
193 write(Ostream& os) const
194 {
196  writeEntry(os, "phi", phiName_);
197  writeEntry(os, "rho", rhoName_);
198  writeEntry(os, "massFluxFraction", massFluxFraction_);
199  writeEntry(os, "value", *this);
200 }
201 
202 
203 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204 
205 namespace Foam
206 {
208  (
211  );
212 
213 }
214 
215 // ************************************************************************* //
Foam::surfaceFields.
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:181
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:663
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const scalarField & magSf() const
Return face area magnitudes.
Definition: fvPatch.C:142
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:175
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:369
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:63
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:243
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.
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.
phi
Definition: correctPhi.H:3
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:351
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
const GeometricField::Patch & lookupPatchField(const word &name) const
Lookup and return the patchField of the named field from the.
Abstract base class for thermophysical transport models (RAS, LES and laminar).
label patchi
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
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:386
const scalarField & deltaCoeffs() const
Return the face - cell distance coeffient.
Definition: fvPatch.C:170
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:216
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:147
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:145
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:357