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 | Copyright (C) 2011-2016 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 
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
31 #include "IOobjectList.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
38 (
39  const fvPatch& p,
41 )
42 :
44  phiName_("phi"),
45  rhoName_("none"),
46  massFluxFraction_(1.0)
47 {
48  refValue() = 0.0;
49  refGrad() = 0.0;
50  valueFraction() = 0.0;
51 }
52 
53 
56 (
57  const fvPatch& p,
59  const dictionary& dict
60 )
61 :
63  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
64  rhoName_(dict.lookupOrDefault<word>("rho", "none")),
65  massFluxFraction_(dict.lookupOrDefault<scalar>("massFluxFraction", 1.0))
66 {
67 
68  refValue() = 1.0;
69  refGrad() = 0.0;
70  valueFraction() = 0.0;
71 
72  if (dict.found("value"))
73  {
75  (
76  Field<scalar>("value", dict, p.size())
77  );
78  }
79  else
80  {
82  }
83 }
84 
87 (
89  const fvPatch& p,
91  const fvPatchFieldMapper& mapper
92 )
93 :
94  mixedFvPatchField<scalar>(ptf, p, iF, mapper),
95  phiName_(ptf.phiName_),
96  rhoName_(ptf.rhoName_),
97  massFluxFraction_(ptf.massFluxFraction_)
98 {}
99 
100 
103 (
105 )
106 :
108  phiName_(tppsf.phiName_),
109  rhoName_(tppsf.rhoName_),
110  massFluxFraction_(tppsf.massFluxFraction_)
111 {}
112 
115 (
118 )
119 :
120  mixedFvPatchField<scalar>(tppsf, iF),
121  phiName_(tppsf.phiName_),
122  rhoName_(tppsf.rhoName_),
123  massFluxFraction_(tppsf.massFluxFraction_)
124 {}
125 
126 
127 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
128 
130 (
131  const fvPatchFieldMapper& m
132 )
133 {
134  scalarField::autoMap(m);
135 }
136 
137 
139 (
140  const fvPatchScalarField& ptf,
141  const labelList& addr
142 )
143 {
145 }
146 
147 
149 {
150  if (this->updated())
151  {
152  return;
153  }
154 
155  const label patchi = patch().index();
156 
158  db().lookupObject
159  <
161  >
162  (
164  (
166  internalField().group()
167  )
168  );
169 
170  const fvsPatchField<scalar>& phip =
171  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
172 
173  const scalarField alphap(turbModel.alphaEff(patchi));
174 
175  refValue() = massFluxFraction_;
176  refGrad() = 0.0;
177 
178  valueFraction() =
179  1.0
180  /
181  (
182  1.0 +
183  alphap*patch().deltaCoeffs()*patch().magSf()/max(mag(phip), SMALL)
184  );
185 
187 
188  if (debug)
189  {
190  scalar phi = gSum(-phip*(*this));
191 
192  Info<< patch().boundaryMesh().mesh().name() << ':'
193  << patch().name() << ':'
194  << this->internalField().name() << " :"
195  << " mass flux[Kg/s]:" << phi
196  << endl;
197  }
198 }
199 
200 
202 write(Ostream& os) const
203 {
205  os.writeKeyword("phi") << phiName_ << token::END_STATEMENT << nl;
206  os.writeKeyword("rho") << rhoName_ << token::END_STATEMENT << nl;
207  os.writeKeyword("massFluxFraction") << massFluxFraction_
208  << token::END_STATEMENT << nl;
209  this->writeEntry("value", os);
210 }
211 
212 
213 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
214 
215 namespace Foam
216 {
218  (
221  );
222 
223 }
224 
225 // ************************************************************************* //
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:339
Foam::surfaceFields.
const char *const group
Group name for atomic constants.
surfaceScalarField & phi
This BC is used for species inlets. The diffusion and advection fluxes are considered to calculate th...
const scalarField & magSf() const
Return face area magnitudes.
Definition: fvPatch.C:136
dictionary dict
Templated abstract base class for LES SGS models.
Definition: LESModel.H:56
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 & deltaCoeffs() const
Return the face - cell distance coeffient.
Definition: fvPatch.C:164
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
const word & name() const
Return name.
Definition: fvPatch.H:149
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:719
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
const fvMesh & mesh() const
Return the mesh reference.
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()
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Type gSum(const FieldField< Field, Type > &f)
static const word propertiesName
Default name of the turbulence properties dictionary.
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.
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:370
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=NULL, const Type *=NULL) const
Lookup and return the patchField of the named field from the.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:262
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:198
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:179
label patchi
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:312
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
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.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
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()
const word & name() const
Return reference to name.
Definition: fvMesh.H:257
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:346
Namespace for OpenFOAM.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:363
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:185