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-2024 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  fixedJumpFvPatchScalarField(p, iF, dict, true),
41  phiName_
42  (
43  cyclicPatch().owner()
44  ? dict.lookupOrDefault<word>("phi", "phi")
45  : word::null
46  ),
47  rhoName_
48  (
49  cyclicPatch().owner()
50  ? dict.lookupOrDefault<word>("rho", "rho")
51  : word::null
52  ),
53  D_
54  (
55  cyclicPatch().owner()
56  ? dict.lookup<scalar>("D", dimless/dimArea)
57  : NaN
58  ),
59  I_
60  (
61  cyclicPatch().owner()
62  ? dict.lookup<scalar>("I", dimless/dimLength)
63  : NaN
64  ),
65  length_
66  (
67  cyclicPatch().owner()
68  ? dict.lookup<scalar>("length", dimLength)
69  : NaN
70  ),
71  relaxation_
72  (
73  cyclicPatch().owner()
74  ? dict.lookupOrDefault<scalar>("relaxation", unitFraction, 1)
75  : NaN
76  ),
77  jump0_(cyclicPatch().owner() ? jump()() : scalarField(p.size()))
78 {}
79 
80 
82 (
84  const fvPatch& p,
86  const fieldMapper& mapper
87 )
88 :
89  fixedJumpFvPatchScalarField(ptf, p, iF, mapper),
90  phiName_(ptf.phiName_),
91  rhoName_(ptf.rhoName_),
92  D_(ptf.D_),
93  I_(ptf.I_),
94  length_(ptf.length_),
95  relaxation_(ptf.relaxation_),
96  jump0_(ptf.jump0_)
97 {}
98 
99 
101 (
104 )
105 :
106  fixedJumpFvPatchScalarField(ptf, iF),
107  phiName_(ptf.phiName_),
108  rhoName_(ptf.rhoName_),
109  D_(ptf.D_),
110  I_(ptf.I_),
111  length_(ptf.length_),
112  relaxation_(ptf.relaxation_),
113  jump0_(ptf.jump0_)
114 {}
115 
116 
117 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
118 
120 {
121  if (updated())
122  {
123  return;
124  }
125 
126  if (cyclicPatch().owner())
127  {
128  const surfaceScalarField& phi =
129  db().lookupObject<surfaceScalarField>(phiName_);
130 
131  const fvsPatchField<scalar>& phip =
132  patch().patchField<surfaceScalarField, scalar>(phi);
133 
134  scalarField Un(phip/patch().magSf());
135 
136  if (phi.dimensions() == dimMassFlux)
137  {
138  Un /= patch().lookupPatchField<volScalarField, scalar>(rhoName_);
139  }
140 
141  const scalarField magUn(mag(Un));
142 
143  const momentumTransportModel& turbModel =
144  db().lookupType<momentumTransportModel>(internalField().group());
145 
146  jumpRef() =
147  -sign(Un)
148  *(
149  D_*turbModel.nu(patch().index())
150  + I_*0.5*magUn
151  )*magUn*length_;
152 
153  if (internalField().dimensions() == dimPressure)
154  {
155  jumpRef() *=
156  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
157  }
158 
159  if (relaxation_ < 1)
160  {
161  jumpRef() += (1 - relaxation_)*(jump0_ - jumpRef());
162  }
163 
164  jump0_ = jumpRef();
165 
166  if (debug)
167  {
168  scalar avePressureJump = gAverage(jumpRef());
169  scalar aveVelocity = gAverage(mag(Un));
170 
171  Info<< patch().boundaryMesh().mesh().name() << ':'
172  << patch().name() << ':'
173  << " Average pressure drop :" << avePressureJump
174  << " Average velocity :" << aveVelocity
175  << endl;
176  }
177  }
178 
179  fixedJumpFvPatchScalarField::updateCoeffs();
180 }
181 
182 
184 {
186 
187  if (cyclicPatch().owner())
188  {
189  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
190  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
191  writeEntry(os, "D", D_);
192  writeEntry(os, "I", I_);
193  writeEntry(os, "length", length_);
194  writeEntryIfDifferent(os, "relaxation", scalar(1), relaxation_);
195  }
196 
197  writeEntry(os, "value", *this);
198 }
199 
200 
201 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
202 
203 namespace Foam
204 {
206  (
209  );
210 }
211 
212 
213 // ************************************************************************* //
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:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
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:83
Abstract base class for turbulence models (RAS, LES and laminar).
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
This boundary condition provides a porous baffle pressure jump condition, using the cyclic condition ...
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
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
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:257
const dimensionSet dimMassFlux
const dimensionSet dimless
SurfaceField< scalar > surfaceScalarField
messageStream Info
const dimensionSet dimLength
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
dimensioned< scalar > mag(const dimensioned< Type > &)
Type gAverage(const FieldField< Field, Type > &f)
const dimensionSet dimArea
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
const unitConversion unitFraction
dictionary dict
volScalarField & p
Foam::surfaceFields.