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 (dict.found("value"))
71  {
73  (
74  Field<scalar>("value", dict, p.size())
75  );
76  }
77  else
78  {
80  }
81  }
82 }
83 
84 
86 (
87  const porousBafflePressureFvPatchField& ptf,
88  const fvPatch& p,
89  const DimensionedField<scalar, volMesh>& iF,
90  const fvPatchFieldMapper& mapper
91 )
92 :
93  fixedJumpFvPatchField<scalar>(ptf, p, iF, mapper),
94  phiName_(ptf.phiName_),
95  rhoName_(ptf.rhoName_),
96  D_(ptf.D_.clone()),
97  I_(ptf.I_.clone()),
98  length_(ptf.length_),
99  uniformJump_(ptf.uniformJump_)
100 {}
101 
102 
104 (
106 )
107 :
109  fixedJumpFvPatchField<scalar>(ptf),
110  phiName_(ptf.phiName_),
111  rhoName_(ptf.rhoName_),
112  D_(ptf.D_.clone()),
113  I_(ptf.I_.clone()),
114  length_(ptf.length_),
115  uniformJump_(ptf.uniformJump_)
116 {}
117 
118 
120 (
123 )
124 :
125  fixedJumpFvPatchField<scalar>(ptf, iF),
126  phiName_(ptf.phiName_),
127  rhoName_(ptf.rhoName_),
128  D_(ptf.D_.clone()),
129  I_(ptf.I_.clone()),
130  length_(ptf.length_),
131  uniformJump_(ptf.uniformJump_)
132 {}
133 
134 
135 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
136 
138 {
139  if (updated())
140  {
141  return;
142  }
143 
144  const auto& phi = db().lookupObject<surfaceScalarField>(phiName_);
145 
146  const fvsPatchField<scalar>& phip =
147  patch().patchField<surfaceScalarField, scalar>(phi);
148 
149  scalarField Un(phip/patch().magSf());
150 
151  if (phi.dimensions() == dimMass/dimTime)
152  {
153  Un /= patch().lookupPatchField<volScalarField, scalar>(rhoName_);
154  }
155 
156  if (uniformJump_)
157  {
158  Un = gAverage(Un);
159  }
160  scalarField magUn(mag(Un));
161 
162  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
163  (
165  (
167  internalField().group()
168  )
169  );
170 
171  const scalar t = db().time().timeOutputValue();
172  const scalar D = D_->value(t);
173  const scalar I = I_->value(t);
174 
175  setJump
176  (
177  -sign(Un)
178  *(
179  D*turbModel.nu(patch().index())
180  + I*0.5*magUn
181  )*magUn*length_
182  );
183 
184  if (internalField().dimensions() == dimPressure)
185  {
186  setJump
187  (
188  jump()*patch().lookupPatchField<volScalarField, scalar>(rhoName_)
189  );
190  }
191 
192  this->relax();
193 
194  if (debug)
195  {
196  scalar avePressureJump = gAverage(jump());
197  scalar aveVelocity = gAverage(Un);
198 
199  Info<< patch().boundaryMesh().mesh().name() << ':'
200  << patch().name() << ':'
201  << " Average pressure drop :" << avePressureJump
202  << " Average velocity :" << aveVelocity
203  << endl;
204  }
205 
207 }
208 
209 
211 {
213  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
214  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
215  D_->writeData(os);
216  I_->writeData(os);
217  os.writeEntry("length", length_);
218  os.writeEntry("uniformJump", uniformJump_);
219 }
220 
221 
222 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
223 
224 namespace Foam
225 {
227  (
229  porousBafflePressureFvPatchField
230  );
231 }
232 
233 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
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:120
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:68
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:312
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
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:84
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:327
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:270
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.