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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30#include "surfaceFields.H"
31#include "turbulenceModel.H"
33
34// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35
37(
38 const fvPatch& p,
40)
41:
42 fixedJumpFvPatchField<scalar>(p, iF),
43 phiName_("phi"),
44 rhoName_("rho"),
45 D_(),
46 I_(),
47 length_(0),
48 uniformJump_(false)
49{}
50
51
53(
54 const fvPatch& p,
56 const dictionary& dict
57)
58:
59 fixedJumpFvPatchField<scalar>(p, iF),
60 phiName_(dict.getOrDefault<word>("phi", "phi")),
61 rhoName_(dict.getOrDefault<word>("rho", "rho")),
62 D_(Function1<scalar>::New("D", dict, &db())),
63 I_(Function1<scalar>::New("I", dict, &db())),
64 length_(dict.get<scalar>("length")),
65 uniformJump_(dict.getOrDefault("uniformJump", false))
66{
68 (
69 Field<scalar>("value", dict, p.size())
70 );
71}
72
73
75(
77 const fvPatch& p,
79 const fvPatchFieldMapper& mapper
80)
81:
82 fixedJumpFvPatchField<scalar>(ptf, p, iF, mapper),
83 phiName_(ptf.phiName_),
84 rhoName_(ptf.rhoName_),
85 D_(ptf.D_.clone()),
86 I_(ptf.I_.clone()),
87 length_(ptf.length_),
88 uniformJump_(ptf.uniformJump_)
89{}
90
91
93(
95)
96:
98 fixedJumpFvPatchField<scalar>(ptf),
99 phiName_(ptf.phiName_),
100 rhoName_(ptf.rhoName_),
101 D_(ptf.D_.clone()),
102 I_(ptf.I_.clone()),
103 length_(ptf.length_),
104 uniformJump_(ptf.uniformJump_)
105{}
106
107
109(
112)
113:
114 fixedJumpFvPatchField<scalar>(ptf, iF),
115 phiName_(ptf.phiName_),
116 rhoName_(ptf.rhoName_),
117 D_(ptf.D_.clone()),
118 I_(ptf.I_.clone()),
119 length_(ptf.length_),
120 uniformJump_(ptf.uniformJump_)
121{}
122
123
124// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
125
127{
128 if (updated())
129 {
130 return;
131 }
132
133 const auto& phi = db().lookupObject<surfaceScalarField>(phiName_);
134
135 const fvsPatchField<scalar>& phip =
136 patch().patchField<surfaceScalarField, scalar>(phi);
137
138 scalarField Un(phip/patch().magSf());
139
140 if (phi.dimensions() == dimMass/dimTime)
141 {
142 Un /= patch().lookupPatchField<volScalarField, scalar>(rhoName_);
143 }
144
145 if (uniformJump_)
146 {
147 Un = gAverage(Un);
148 }
149 scalarField magUn(mag(Un));
150
151 const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
152 (
154 (
156 internalField().group()
157 )
158 );
159
160 const scalar t = db().time().timeOutputValue();
161 const scalar D = D_->value(t);
162 const scalar I = I_->value(t);
163
164 setJump
165 (
166 -sign(Un)
167 *(
168 D*turbModel.nu(patch().index())
169 + I*0.5*magUn
170 )*magUn*length_
171 );
172
173 if (internalField().dimensions() == dimPressure)
174 {
175 setJump
176 (
177 jump()*patch().lookupPatchField<volScalarField, scalar>(rhoName_)
178 );
179 }
180
181 if (debug)
182 {
183 scalar avePressureJump = gAverage(jump());
184 scalar aveVelocity = gAverage(Un);
185
186 Info<< patch().boundaryMesh().mesh().name() << ':'
187 << patch().name() << ':'
188 << " Average pressure drop :" << avePressureJump
189 << " Average velocity :" << aveVelocity
190 << endl;
191 }
192
194}
195
196
198{
200 os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
201 os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
202 D_->writeData(os);
203 I_->writeData(os);
204 os.writeEntry("length", length_);
205 os.writeEntry("uniformJump", uniformJump_);
206}
207
208
209// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
210
211namespace Foam
212{
214 (
217 );
218}
219
220// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
surfaceScalarField & phi
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
Generic templated field type.
Definition: Field.H:82
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:96
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:251
scalar timeOutputValue() const
Return current time value.
Definition: TimeStateI.H:31
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Abstract base class for cyclic coupled interfaces.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const Type & value() const
Return const reference to value.
This boundary condition provides a jump condition, using the cyclic condition as a base.
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:79
This boundary condition provides a jump condition, using the cyclic condition as a base.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Abstract base class for turbulence models (RAS, LES and laminar).
static const word propertiesName
Default name of the turbulence properties dictionary.
const Time & time() const
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
Namespace for OpenFOAM.
const dimensionSet dimPressure
dimensionedScalar sign(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
messageStream Info
Information stream (stdout output on master, null elsewhere)
static const Identity< scalar > I
Definition: Identity.H:94
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Type gAverage(const FieldField< Field, Type > &f)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
dictionary dict
const dimensionedScalar & D
Foam::surfaceFields.