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 -------------------------------------------------------------------------------
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 )
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 
211 namespace Foam
212 {
214  (
216  porousBafflePressureFvPatchField
217  );
218 }
219 
220 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:47
Foam::fvPatchScalarField
fvPatchField< scalar > fvPatchScalarField
Definition: fvPatchFieldsFwd.H:33
Foam::TimeState::timeOutputValue
scalar timeOutputValue() const
Definition: TimeStateI.H:24
Foam::dimPressure
const dimensionSet dimPressure
Foam::fvPatchField< scalar >::internalField
const DimensionedField< Type, volMesh > & internalField() const
Definition: fvPatchField.H:359
Foam::porousBafflePressureFvPatchField::write
virtual void write(Ostream &) const
Definition: porousBafflePressureFvPatchField.C:190
Foam::Ostream::writeEntryIfDifferent
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Definition: Ostream.H:244
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:46
Foam::turbulenceModel::nu
virtual tmp< volScalarField > nu() const =0
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::fixedJumpFvPatchField< scalar >::setJump
virtual void setJump(const Field< scalar > &jump)
Definition: fixedJumpFvPatchField.C:140
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::constant::atomic::group
constexpr const char *const group
Definition: atomicConstants.H:46
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:597
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:64
Foam::turbulenceModel::propertiesName
static const word propertiesName
Definition: turbulenceModel.H:96
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
surfaceFields.H
Foam::surfaceFields.
Foam::dimensioned::value
const Type & value() const
Definition: dimensionedType.C:427
Foam::fvPatch::name
virtual const word & name() const
Definition: fvPatch.H:163
Foam::fvBoundaryMesh::mesh
const fvMesh & mesh() const noexcept
Definition: fvBoundaryMesh.H:99
Foam::sign
dimensionedScalar sign(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:159
Foam::fvPatchField< scalar >::updated
bool updated() const
Definition: fvPatchField.H:383
Foam::Function1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: propellerInfo.H:287
Foam::fvPatch::boundaryMesh
const fvBoundaryMesh & boundaryMesh() const
Definition: fvPatch.H:199
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::Field< scalar >
Foam::compressible::turbulenceModel
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Definition: turbulentFluidThermoModel.H:59
Foam::Info
messageStream Info
Foam::fixedJumpFvPatchField< scalar >
Foam::fixedJumpFvPatchField::write
virtual void write(Ostream &) const
Definition: fixedJumpFvPatchField.C:247
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:162
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:53
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:427
Foam::fvPatch::lookupPatchField
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const Type *=nullptr) const
Definition: fvPatchFvMeshTemplates.C:27
Foam::porousBafflePressureFvPatchField::porousBafflePressureFvPatchField
porousBafflePressureFvPatchField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Definition: porousBafflePressureFvPatchField.C:30
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:59
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:119
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Definition: atmBoundaryLayer.C:26
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Foam::fvPatchField< Type >::updateCoeffs
virtual void updateCoeffs()
Definition: fvPatchField.C:314
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Definition: DimensionedFieldReuseFunctions.H:100
Foam::porousBafflePressureFvPatchField::updateCoeffs
virtual void updateCoeffs()
Definition: porousBafflePressureFvPatchField.C:119
Foam::cyclicLduInterfaceField
Abstract base class for cyclic coupled interfaces.
Definition: cyclicLduInterfaceField.H:48
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:50
Foam::fvPatchField< scalar >::db
const objectRegistry & db() const
Definition: fvPatchField.C:199
Foam::foamVersion::patch
const std::string patch
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
D
const dimensionedScalar & D
Definition: solveBulkSurfactant.H:4
Foam::fvPatch::patchField
const GeometricField::Patch & patchField(const GeometricField &) const
Definition: fvPatchTemplates.C:75
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Definition: Ostream.H:232
Foam::PtrListOps::get
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
porousBafflePressureFvPatchField.H
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:41
Foam::fvPatchField< scalar >::patch
const fvPatch & patch() const
Definition: fvPatchField.H:353
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:52
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:49
turbulenceModel.H
Foam::objectRegistry::time
const Time & time() const noexcept
Definition: objectRegistry.H:174
Foam::I
static const Identity< scalar > I
Definition: Identity.H:89
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:50
Foam::fixedJumpFvPatchField< scalar >::jump
virtual tmp< Field< scalar > > jump() const
Definition: fixedJumpFvPatchField.C:160
Foam::fvMesh::name
const word & name() const
Definition: fvMesh.H:296