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-2021 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 "momentumTransportModel.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_(dict.lookup<scalar>("D")),
59  I_(dict.lookup<scalar>("I")),
60  length_(dict.lookup<scalar>("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 (
90 )
91 :
93  phiName_(ptf.phiName_),
94  rhoName_(ptf.rhoName_),
95  D_(ptf.D_),
96  I_(ptf.I_),
97  length_(ptf.length_)
98 {}
99 
100 
101 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
102 
104 {
105  if (updated())
106  {
107  return;
108  }
109 
110  const surfaceScalarField& phi =
111  db().lookupObject<surfaceScalarField>(phiName_);
112 
113  const fvsPatchField<scalar>& phip =
115 
116  scalarField Un(phip/patch().magSf());
117 
118  if (phi.dimensions() == dimMassFlux)
119  {
120  Un /= patch().lookupPatchField<volScalarField, scalar>(rhoName_);
121  }
122 
123  scalarField magUn(mag(Un));
124 
125  const momentumTransportModel& turbModel =
127  (
129  (
130  momentumTransportModel::typeName,
131  internalField().group()
132  )
133  );
134 
135  jump_ =
136  -sign(Un)
137  *(
138  D_*turbModel.nu(patch().index())
139  + I_*0.5*magUn
140  )*magUn*length_;
141 
142  if (internalField().dimensions() == dimPressure)
143  {
144  jump_ *= patch().lookupPatchField<volScalarField, scalar>(rhoName_);
145  }
146 
147  if (debug)
148  {
149  scalar avePressureJump = gAverage(jump_);
150  scalar aveVelocity = gAverage(mag(Un));
151 
152  Info<< patch().boundaryMesh().mesh().name() << ':'
153  << patch().name() << ':'
154  << " Average pressure drop :" << avePressureJump
155  << " Average velocity :" << aveVelocity
156  << endl;
157  }
158 
160 }
161 
162 
164 {
166  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
167  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
168  writeEntry(os, "D", D_);
169  writeEntry(os, "I", I_);
170  writeEntry(os, "length", length_);
171 }
172 
173 
174 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
175 
176 namespace Foam
177 {
179  (
182  );
183 }
184 
185 // ************************************************************************* //
Foam::surfaceFields.
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:181
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.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:369
const dimensionSet dimPressure
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
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.
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:58
const dimensionSet & dimensions() const
Return dimensions.
const fvMesh & mesh() const
Return the mesh reference.
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
phi
Definition: correctPhi.H:3
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:351
virtual void write(Ostream &) const
Write.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const GeometricField::Patch & lookupPatchField(const word &name) const
Lookup and return the patchField of the named field from the.
Abstract base class for turbulence models (RAS, LES and laminar).
This boundary condition provides a jump condition, using the cyclic condition as a base...
const dimensionSet dimMassFlux
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:386
Type gAverage(const FieldField< Field, Type > &f)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:216
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:147
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:145
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
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:357