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-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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  const bool valueRequired
58 )
59 :
60  fixedJumpFvPatchField<scalar>(p, iF, dict, false),
61  phiName_(dict.getOrDefault<word>("phi", "phi")),
62  rhoName_(dict.getOrDefault<word>("rho", "rho")),
63  D_(Function1<scalar>::New("D", dict, &db())),
64  I_(Function1<scalar>::New("I", dict, &db())),
65  length_(dict.get<scalar>("length")),
66  uniformJump_(dict.getOrDefault("uniformJump", false))
67 {
68  if (valueRequired)
69  {
70  if (!this->readValueEntry(dict))
71  {
73  }
74  }
75 }
76 
77 
79 (
80  const porousBafflePressureFvPatchField& ptf,
81  const fvPatch& p,
82  const DimensionedField<scalar, volMesh>& iF,
83  const fvPatchFieldMapper& mapper
84 )
85 :
86  fixedJumpFvPatchField<scalar>(ptf, p, iF, mapper),
87  phiName_(ptf.phiName_),
88  rhoName_(ptf.rhoName_),
89  D_(ptf.D_.clone()),
90  I_(ptf.I_.clone()),
91  length_(ptf.length_),
92  uniformJump_(ptf.uniformJump_)
93 {}
94 
95 
97 (
99 )
100 :
102  fixedJumpFvPatchField<scalar>(ptf),
103  phiName_(ptf.phiName_),
104  rhoName_(ptf.rhoName_),
105  D_(ptf.D_.clone()),
106  I_(ptf.I_.clone()),
107  length_(ptf.length_),
108  uniformJump_(ptf.uniformJump_)
109 {}
110 
111 
113 (
116 )
117 :
118  fixedJumpFvPatchField<scalar>(ptf, iF),
119  phiName_(ptf.phiName_),
120  rhoName_(ptf.rhoName_),
121  D_(ptf.D_.clone()),
122  I_(ptf.I_.clone()),
123  length_(ptf.length_),
124  uniformJump_(ptf.uniformJump_)
125 {}
126 
127 
128 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
129 
131 {
132  if (updated())
133  {
134  return;
135  }
136 
137  const auto& phip = patch().lookupPatchField<surfaceScalarField>(phiName_);
138 
139  scalarField Un(phip/patch().magSf());
140 
141  if (phip.internalField().dimensions() == dimMass/dimTime)
142  {
143  Un /= patch().lookupPatchField<volScalarField>(rhoName_);
144  }
145 
146  if (uniformJump_)
147  {
148  Un = gAverage(Un);
149  }
150  scalarField magUn(mag(Un));
151 
152  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
153  (
155  (
157  internalField().group()
158  )
159  );
160 
161  const scalar t = db().time().timeOutputValue();
162  const scalar D = D_->value(t);
163  const scalar I = I_->value(t);
164 
165  setJump
166  (
167  -sign(Un)
168  *(
169  D*turbModel.nu(patch().index())
170  + I*0.5*magUn
171  )*magUn*length_
172  );
173 
174  if (internalField().dimensions() == dimPressure)
175  {
176  setJump
177  (
178  jump()*patch().lookupPatchField<volScalarField>(rhoName_)
179  );
180  }
181 
182  this->relax();
183 
184  if (debug)
185  {
186  scalar avePressureJump = gAverage(jump());
187  scalar aveVelocity = gAverage(Un);
188 
189  Info<< patch().boundaryMesh().mesh().name() << ':'
190  << patch().name() << ':'
191  << " Average pressure drop :" << avePressureJump
192  << " Average velocity :" << aveVelocity
193  << endl;
194  }
195 
197 }
198 
199 
201 {
203  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
204  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
205  D_->writeData(os);
206  I_->writeData(os);
207  os.writeEntry("length", length_);
208  os.writeEntry("uniformJump", uniformJump_);
209 }
210 
211 
212 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
213 
214 namespace Foam
215 {
217  (
219  porousBafflePressureFvPatchField
220  );
221 }
222 
223 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the "value" entry into *this.
Definition: fvPatchField.C:32
Foam::surfaceFields.
dimensionedScalar sign(const dimensionedScalar &ds)
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
const Type & value() const noexcept
Return const reference to value.
dictionary dict
"blocking" : (MPI_Bsend, MPI_Recv)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
UEqn relax()
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Macros for easy insertion into run-time selection tables.
constexpr const char *const group
Group name for atomic constants.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
static const Identity< scalar > I
Definition: Identity.H:100
fvPatchField< scalar > fvPatchScalarField
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual void evaluate(const Pstream::commsTypes commsType)
Evaluate the patch field.
const dimensionSet dimPressure
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:336
virtual void write(Ostream &) const
Write.
int debug
Static debugging option.
OBJstream os(runTime.globalPath()/outputName)
This boundary condition provides a jump condition, using the cyclic condition as a base...
This boundary condition provides a jump condition, using the cyclic condition as a base...
Type gAverage(const FieldField< Field, Type > &f)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:309
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
porousBafflePressureFvPatchField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
messageStream Info
Information stream (stdout output on master, null elsewhere)
const dimensionedScalar & D
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Abstract base class for cyclic coupled interfaces.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.