porousBafflePressureFvPatchField.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-2015 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 :
93  phiName_(ptf.phiName_),
94  rhoName_(ptf.rhoName_),
95  D_(ptf.D_),
96  I_(ptf.I_),
97  length_(ptf.length_)
98 {}
99 
100 
102 (
105 )
106 :
108  phiName_(ptf.phiName_),
109  rhoName_(ptf.rhoName_),
110  D_(ptf.D_),
111  I_(ptf.I_),
112  length_(ptf.length_)
113 {}
114 
115 
116 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
117 
119 {
120  if (updated())
121  {
122  return;
123  }
124 
125  const surfaceScalarField& phi =
126  db().lookupObject<surfaceScalarField>(phiName_);
127 
128  const fvsPatchField<scalar>& phip =
130 
131  scalarField Un(phip/patch().magSf());
132 
134  {
135  Un /= patch().lookupPatchField<volScalarField, scalar>(rhoName_);
136  }
137 
138  scalarField magUn(mag(Un));
139 
140  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
141  (
143  (
146  )
147  );
148 
149  jump_ =
150  -sign(Un)
151  *(
152  D_*turbModel.nu(patch().index())
153  + I_*0.5*magUn
154  )*magUn*length_;
155 
156  if (dimensionedInternalField().dimensions() == dimPressure)
157  {
158  jump_ *= patch().lookupPatchField<volScalarField, scalar>(rhoName_);
159  }
160 
161  if (debug)
162  {
163  scalar avePressureJump = gAverage(jump_);
164  scalar aveVelocity = gAverage(mag(Un));
165 
166  Info<< patch().boundaryMesh().mesh().name() << ':'
167  << patch().name() << ':'
168  << " Average pressure drop :" << avePressureJump
169  << " Average velocity :" << aveVelocity
170  << endl;
171  }
172 
174 }
175 
176 
178 {
180  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
181  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
182  os.writeKeyword("D") << D_ << token::END_STATEMENT << nl;
183  os.writeKeyword("I") << I_ << token::END_STATEMENT << nl;
184  os.writeKeyword("length") << length_ << token::END_STATEMENT << nl;
185 }
186 
187 
188 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
189 
190 namespace Foam
191 {
193  (
196  );
197 }
198 
199 // ************************************************************************* //
const dimensionSet dimPressure
Foam::surfaceFields.
const word & name() const
Return reference to name.
Definition: fvMesh.H:257
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
#define readScalar
Definition: doubleScalar.C:38
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:301
const char *const group
Group name for atomic constants.
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
A class for handling words, derived from string.
Definition: word.H:59
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
const GeometricField::PatchFieldType & patchField(const GeometricField &) const
Return the corresponding patchField of the named field.
Foam::fvPatchFieldMapper.
dimensionedScalar sign(const dimensionedScalar &ds)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
messageStream Info
const dimensionSet dimDensity
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:348
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Namespace for OpenFOAM.
const word & name() const
Return name.
Definition: fvPatch.H:149
const GeometricField::PatchFieldType & lookupPatchField(const word &name, const GeometricField *=NULL, const Type *=NULL) const
Lookup and return the patchField of the named field from the.
Abstract base class for cyclic coupled interfaces.
static word groupName(Name name, const word &group)
dictionary dict
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:185
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
volScalarField & p
Definition: createFields.H:51
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Macros for easy insertion into run-time selection tables.
const dimensionSet & dimensions() const
Return dimensions.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
surfaceScalarField & phi
This boundary condition provides a jump condition, using the cyclic condition as a base...
Abstract base class for turbulence models (RAS, LES and laminar).
const fvMesh & mesh() const
Return the mesh reference.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
Type gAverage(const FieldField< Field, Type > &f)
porousBafflePressureFvPatchField(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,.
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:188
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
static const word propertiesName
Default name of the turbulence properties dictionary.
const DimensionedField< Type, volMesh > & dimensionedInternalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:307
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:300
const dimensionSet dimVelocity
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65