porousBafflePressureFvPatchField.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-2019 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"
28 #include "turbulenceModel.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
34 (
35  const fvPatch& p,
37 )
38 :
40  phiName_("phi"),
41  rhoName_("rho"),
42  D_(0),
43  I_(0),
44  length_(0)
45 {}
46 
47 
49 (
50  const fvPatch& p,
52  const dictionary& dict
53 )
54 :
56  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
57  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
58  D_(readScalar(dict.lookup("D"))),
59  I_(readScalar(dict.lookup("I"))),
60  length_(readScalar(dict.lookup("length")))
61 {
63  (
64  Field<scalar>("value", dict, p.size())
65  );
66 }
67 
68 
70 (
72  const fvPatch& p,
74  const fvPatchFieldMapper& mapper
75 )
76 :
77  fixedJumpFvPatchField<scalar>(ptf, p, iF, mapper),
78  phiName_(ptf.phiName_),
79  rhoName_(ptf.rhoName_),
80  D_(ptf.D_),
81  I_(ptf.I_),
82  length_(ptf.length_)
83 {}
84 
85 
87 (
89 )
90 :
92  phiName_(ptf.phiName_),
93  rhoName_(ptf.rhoName_),
94  D_(ptf.D_),
95  I_(ptf.I_),
96  length_(ptf.length_)
97 {}
98 
99 
101 (
104 )
105 :
107  phiName_(ptf.phiName_),
108  rhoName_(ptf.rhoName_),
109  D_(ptf.D_),
110  I_(ptf.I_),
111  length_(ptf.length_)
112 {}
113 
114 
115 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116 
118 {
119  if (updated())
120  {
121  return;
122  }
123 
124  const surfaceScalarField& phi =
125  db().lookupObject<surfaceScalarField>(phiName_);
126 
127  const fvsPatchField<scalar>& phip =
129 
130  scalarField Un(phip/patch().magSf());
131 
133  {
134  Un /= patch().lookupPatchField<volScalarField, scalar>(rhoName_);
135  }
136 
137  scalarField magUn(mag(Un));
138 
139  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
140  (
142  (
144  internalField().group()
145  )
146  );
147 
148  jump_ =
149  -sign(Un)
150  *(
151  D_*turbModel.nu(patch().index())
152  + I_*0.5*magUn
153  )*magUn*length_;
154 
155  if (internalField().dimensions() == dimPressure)
156  {
157  jump_ *= patch().lookupPatchField<volScalarField, scalar>(rhoName_);
158  }
159 
160  if (debug)
161  {
162  scalar avePressureJump = gAverage(jump_);
163  scalar aveVelocity = gAverage(mag(Un));
164 
165  Info<< patch().boundaryMesh().mesh().name() << ':'
166  << patch().name() << ':'
167  << " Average pressure drop :" << avePressureJump
168  << " Average velocity :" << aveVelocity
169  << endl;
170  }
171 
173 }
174 
175 
177 {
179  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
180  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
181  writeEntry(os, "D", D_);
182  writeEntry(os, "I", I_);
183  writeEntry(os, "length", length_);
184 }
185 
186 
187 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 
189 namespace Foam
190 {
192  (
195  );
196 }
197 
198 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
Foam::surfaceFields.
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:179
const char *const group
Group name for atomic constants.
dimensionedScalar sign(const dimensionedScalar &ds)
dictionary dict
const GeometricField::Patch & patchField(const GeometricField &) const
Return the corresponding patchField of the named field.
surfaceScalarField & phi
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:356
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Abstract base class for turbulence models (RAS, LES and laminar).
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
const dimensionSet & dimensions() const
Return dimensions.
const fvMesh & mesh() const
Return the mesh reference.
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)
Foam::fvPatchFieldMapper.
const dimensionSet dimPressure
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:325
virtual void write(Ostream &) const
Write.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const dimensionSet dimDensity
This boundary condition provides a jump condition, using the cyclic condition as a base...
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:253
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const Type *=nullptr) const
Lookup and return the patchField of the named field from the.
Type gAverage(const FieldField< Field, Type > &f)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:229
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:167
porousBafflePressureFvPatchField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual const word & name() const
Return name.
Definition: fvPatch.H:143
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
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
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:583
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:332
const dimensionSet dimVelocity