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_(),
43  I_(),
44  length_(0),
45  uniformJump_(false)
46 {}
47 
48 
50 (
51  const fvPatch& p,
53  const dictionary& dict
54 )
55 :
57  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
58  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
59  D_(DataEntry<scalar>::New("D", dict)),
60  I_(DataEntry<scalar>::New("I", dict)),
61  length_(readScalar(dict.lookup("length"))),
62  uniformJump_(dict.lookupOrDefault<bool>("uniformJump", false))
63 {
65  (
66  Field<scalar>("value", dict, p.size())
67  );
68 }
69 
70 
72 (
74  const fvPatch& p,
76  const fvPatchFieldMapper& mapper
77 )
78 :
79  fixedJumpFvPatchField<scalar>(ptf, p, iF, mapper),
80  phiName_(ptf.phiName_),
81  rhoName_(ptf.rhoName_),
82  D_(ptf.D_, false),
83  I_(ptf.I_, false),
84  length_(ptf.length_),
85  uniformJump_(ptf.uniformJump_)
86 {}
87 
88 
90 (
92 )
93 :
96  phiName_(ptf.phiName_),
97  rhoName_(ptf.rhoName_),
98  D_(ptf.D_, false),
99  I_(ptf.I_, false),
100  length_(ptf.length_),
101  uniformJump_(ptf.uniformJump_)
102 {}
103 
104 
106 (
109 )
110 :
112  phiName_(ptf.phiName_),
113  rhoName_(ptf.rhoName_),
114  D_(ptf.D_, false),
115  I_(ptf.I_, false),
116  length_(ptf.length_),
117  uniformJump_(ptf.uniformJump_)
118 {}
119 
120 
121 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
122 
124 {
125  if (updated())
126  {
127  return;
128  }
129 
130  const surfaceScalarField& phi =
132 
133  const fvsPatchField<scalar>& phip =
135 
136  scalarField Un(phip/patch().magSf());
137 
138  if (phi.dimensions() == dimMass/dimTime)
139  {
140  Un /= patch().lookupPatchField<volScalarField, scalar>(rhoName_);
141  }
142 
143  if (uniformJump_)
144  {
145  Un = gAverage(Un);
146  }
147  scalarField magUn(mag(Un));
148 
149  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
150  (
152  (
155  )
156  );
157 
158  const scalar t = db().time().timeOutputValue();
159  const scalar D = D_->value(t);
160  const scalar I = I_->value(t);
161 
162  jump_ =
163  -sign(Un)
164  *(
165  D*turbModel.nu(patch().index())
166  + I*0.5*magUn
167  )*magUn*length_;
168 
169  if (dimensionedInternalField().dimensions() == dimPressure)
170  {
172  }
173 
174  if (debug)
175  {
176  scalar avePressureJump = gAverage(jump_);
177  scalar aveVelocity = gAverage(Un);
178 
179  Info<< patch().boundaryMesh().mesh().name() << ':'
180  << patch().name() << ':'
181  << " Average pressure drop :" << avePressureJump
182  << " Average velocity :" << aveVelocity
183  << endl;
184  }
185 
187 }
188 
189 
191 {
193  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
194  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
195  D_->writeData(os);
196  I_->writeData(os);
197  os.writeKeyword("length") << length_ << token::END_STATEMENT << nl;
198  os.writeKeyword("uniformJump") << uniformJump_
199  << token::END_STATEMENT << nl;
200 }
201 
202 
203 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204 
205 namespace Foam
206 {
208  (
211  );
212 }
213 
214 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Foam::TimeState::timeOutputValue
scalar timeOutputValue() const
Return current time value.
Definition: TimeState.C:67
Foam::dimPressure
const dimensionSet dimPressure
Foam::porousBafflePressureFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: porousBafflePressureFvPatchField.C:190
Foam::token::END_STATEMENT
@ END_STATEMENT
Definition: token.H:99
Foam::turbulenceModel::nu
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
Foam::porousBafflePressureFvPatchField::I_
autoPtr< DataEntry< scalar > > I_
Inertia pressure lost coefficient.
Definition: porousBafflePressureFvPatchField.H:180
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:571
Foam::fvPatchField< scalar >::dimensionedInternalField
const DimensionedField< Type, volMesh > & dimensionedInternalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:307
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::objectRegistry::time
const Time & time() const
Return time.
Definition: objectRegistry.H:117
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:97
Foam::porousBafflePressureFvPatchField::phiName_
const word phiName_
Name of flux field (default = phi)
Definition: porousBafflePressureFvPatchField.H:171
Foam::constant::atomic::group
const char *const group
Group name for atomic constants.
Definition: atomicConstants.C:40
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceFields.H
Foam::surfaceFields.
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::fvBoundaryMesh::mesh
const fvMesh & mesh() const
Return the mesh reference.
Definition: fvBoundaryMesh.H:99
Foam::sign
dimensionedScalar sign(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:179
Foam::fvPatchField< scalar >::updated
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:348
Foam::fvPatch::boundaryMesh
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:185
Foam::fvPatch::name
const word & name() const
Return name.
Definition: fvPatch.H:149
Foam::porousBafflePressureFvPatchField::rhoName_
const word rhoName_
Name of density field (default = rho)
Definition: porousBafflePressureFvPatchField.H:174
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::porousBafflePressureFvPatchField::length_
scalar length_
Porous media length.
Definition: porousBafflePressureFvPatchField.H:183
Foam::Field< scalar >
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::fixedJumpFvPatchField< scalar >
Foam::fixedJumpFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fixedJumpFvPatchField.C:155
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::porousBafflePressureFvPatchField
This boundary condition provides a jump condition, using the cyclic condition as a base.
Definition: porousBafflePressureFvPatchField.H:164
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Foam::I
static const sphericalTensor I(1)
Foam::porousBafflePressureFvPatchField::porousBafflePressureFvPatchField
porousBafflePressureFvPatchField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: porousBafflePressureFvPatchField.C:34
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:60
Foam::fvPatch::patchField
const GeometricField::PatchFieldType & patchField(const GeometricField &) const
Return the corresponding patchField of the named field.
Definition: fvPatchTemplates.C:70
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Foam::makePatchTypeField
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
Foam::fvPatchField< Type >::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:296
readScalar
#define readScalar
Definition: doubleScalar.C:38
Foam::porousBafflePressureFvPatchField::D_
autoPtr< DataEntry< scalar > > D_
Darcy pressure loss coefficient.
Definition: porousBafflePressureFvPatchField.H:177
Foam::porousBafflePressureFvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: porousBafflePressureFvPatchField.C:123
Foam::cyclicLduInterfaceField
Abstract base class for cyclic coupled interfaces.
Definition: cyclicLduInterfaceField.H:49
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
Foam::fvPatchField< scalar >::db
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:181
Foam::Ostream::writeKeyword
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
Foam::fvPatch::lookupPatchField
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.
Definition: fvPatchFvMeshTemplates.C:32
Foam::porousBafflePressureFvPatchField::uniformJump_
bool uniformJump_
Aplies uniform pressure drop.
Definition: porousBafflePressureFvPatchField.H:186
porousBafflePressureFvPatchField.H
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::fvPatchField< scalar >::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:300
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Definition: objectRegistryTemplates.C:165
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::DataEntry
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: DataEntry.H:52
turbulenceModel.H
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51
Foam::fvMesh::name
const word & name() const
Return reference to name.
Definition: fvMesh.H:257
Foam::fixedJumpFvPatchField< scalar >::jump_
Field< scalar > jump_
"jump" field
Definition: fixedJumpFvPatchField.H:107