Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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
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)
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(
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 = -gWeightedSum(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_);
204
205
206// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
207
208namespace Foam
209{
211 (
213 totalFlowRateAdvectiveDiffusiveFvPatchScalarField
214 );
215
216}
217
218// ************************************************************************* //
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...
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const objectRegistry & db() const
The associated objectRegistry.
const fvPatch & patch() const noexcept
Return the patch.
bool updated() const noexcept
True if the boundary condition has already been updated.
A FieldMapper for finite-volume patch fields.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
const DimensionedField< scalar, volMesh > & internalField() const noexcept
virtual void operator=(const UList< Type > &)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
This boundary condition provides a base class for 'mixed' type boundary conditions,...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
mixedFvPatchField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
virtual Field< scalar > & refGrad()
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual Field< scalar > & refValue()
virtual scalarField & valueFraction()
This BC is used for species inlets. The diffusion and advection fluxes are considered to calculate th...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
totalFlowRateAdvectiveDiffusiveFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< label > labelList
A List of labels.
Definition List.H:62
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Type gWeightedSum(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted sum (integral) of a field, using the mag() of the weights.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fvPatchField< scalar > fvPatchScalarField
dictionary dict
Foam::surfaceFields.