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-2023 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  const dictionary& dict
38 )
39 :
40  fixedJumpFvPatchField<scalar>(p, iF),
41  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
42  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
43  D_(dict.lookup<scalar>("D")),
44  I_(dict.lookup<scalar>("I")),
45  length_(dict.lookup<scalar>("length")),
46  relaxation_(dict.lookupOrDefault<scalar>("relaxation", 1)),
47  jump0_(jump_)
48 {
50  (
51  Field<scalar>("value", dict, p.size())
52  );
53 }
54 
55 
57 (
59  const fvPatch& p,
61  const fvPatchFieldMapper& mapper
62 )
63 :
64  fixedJumpFvPatchField<scalar>(ptf, p, iF, mapper),
65  phiName_(ptf.phiName_),
66  rhoName_(ptf.rhoName_),
67  D_(ptf.D_),
68  I_(ptf.I_),
69  length_(ptf.length_),
70  relaxation_(ptf.relaxation_),
71  jump0_(ptf.jump_)
72 {}
73 
74 
76 (
79 )
80 :
81  fixedJumpFvPatchField<scalar>(ptf, iF),
82  phiName_(ptf.phiName_),
83  rhoName_(ptf.rhoName_),
84  D_(ptf.D_),
85  I_(ptf.I_),
86  length_(ptf.length_),
87  relaxation_(ptf.relaxation_),
88  jump0_(ptf.jump_)
89 {}
90 
91 
92 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
93 
95 {
96  if (updated())
97  {
98  return;
99  }
100 
101  const surfaceScalarField& phi =
102  db().lookupObject<surfaceScalarField>(phiName_);
103 
104  const fvsPatchField<scalar>& phip =
105  patch().patchField<surfaceScalarField, scalar>(phi);
106 
107  scalarField Un(phip/patch().magSf());
108 
109  if (phi.dimensions() == dimMassFlux)
110  {
111  Un /= patch().lookupPatchField<volScalarField, scalar>(rhoName_);
112  }
113 
114  const scalarField magUn(mag(Un));
115 
116  const momentumTransportModel& turbModel =
117  db().lookupType<momentumTransportModel>(internalField().group());
118 
119  jump_ =
120  -sign(Un)
121  *(
122  D_*turbModel.nu(patch().index())
123  + I_*0.5*magUn
124  )*magUn*length_;
125 
126  if (internalField().dimensions() == dimPressure)
127  {
128  jump_ *= patch().lookupPatchField<volScalarField, scalar>(rhoName_);
129  }
130 
131  if (relaxation_ < 1)
132  {
133  jump_ += (1 - relaxation_)*(jump0_ - jump_);
134  }
135 
136  jump0_ = jump_;
137 
138  if (debug)
139  {
140  scalar avePressureJump = gAverage(jump_);
141  scalar aveVelocity = gAverage(mag(Un));
142 
143  Info<< patch().boundaryMesh().mesh().name() << ':'
144  << patch().name() << ':'
145  << " Average pressure drop :" << avePressureJump
146  << " Average velocity :" << aveVelocity
147  << endl;
148  }
149 
151 }
152 
153 
155 {
157  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
158  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
159  writeEntry(os, "D", D_);
160  writeEntry(os, "I", I_);
161  writeEntry(os, "length", length_);
162  writeEntryIfDifferent(os, "relaxation", scalar(1), relaxation_);
163 }
164 
165 
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167 
168 namespace Foam
169 {
171  (
174  );
175 }
176 
177 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
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 GeometricField class.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
This boundary condition provides a jump condition, using the cyclic condition as a base.
virtual void write(Ostream &) const
Write.
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:204
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:82
Abstract base class for turbulence models (RAS, LES and laminar).
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
This boundary condition provides a jump condition, using the cyclic condition as a base.
porousBafflePressureFvPatchField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
void writeEntryIfDifferent(Ostream &os, const word &entryName, const EntryType &value1, const EntryType &value2)
Helper function to write the keyword and entry only if the.
dimensionedScalar sign(const dimensionedScalar &ds)
const dimensionSet dimPressure
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const dimensionSet dimMassFlux
SurfaceField< scalar > surfaceScalarField
messageStream Info
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
dimensioned< scalar > mag(const dimensioned< Type > &)
Type gAverage(const FieldField< Field, Type > &f)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p
Foam::surfaceFields.