omegaWallFunctionFvPatchScalarField.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 
28 #include "momentumTransportModel.H"
29 #include "fvMatrix.H"
31 
32 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
35 Foam::omegaWallFunctionFvPatchScalarField::calculate
36 (
37  const momentumTransportModel& mtm
38 ) const
39 {
40  const label patchi = patch().index();
41 
42  const nutWallFunctionFvPatchScalarField& nutw =
44 
45  const scalarField& y = mtm.yb()[patchi];
46 
47  const tmp<scalarField> tnuw = mtm.nu(patchi);
48  const scalarField& nuw = tnuw();
49 
50  const tmp<volScalarField> tk = mtm.k();
51  const volScalarField& k = tk();
52 
53  const fvPatchVectorField& Uw = mtm.U().boundaryField()[patchi];
54 
55  const scalarField magGradUw(mag(Uw.snGrad()));
56 
57  const VolInternalField<scalar>& Gvf =
58  db().lookupObject<VolInternalField<scalar>>(mtm.GName());
59 
60  const scalar Cmu25 = pow025(nutw.Cmu());
61  const scalar Cmu5 = sqrt(nutw.Cmu());
62 
63  // Initialise the returned field pair
64  Pair<tmp<scalarField>> GandOmega
65  (
66  tmp<scalarField>(new scalarField(patch().size())),
67  tmp<scalarField>(new scalarField(patch().size()))
68  );
69  scalarField& G = GandOmega.first().ref();
70  scalarField& omega = GandOmega.second().ref();
71 
72  // Calculate G and omega
73  forAll(nutw, facei)
74  {
75  const label celli = patch().faceCells()[facei];
76 
77  const scalar Rey = y[facei]*sqrt(k[celli])/nuw[facei];
78  const scalar yPlus = Cmu25*Rey;
79  const scalar uPlus = (1/nutw.kappa())*log(nutw.E()*yPlus);
80 
81  if (blended_)
82  {
83  const scalar lamFrac = exp(-Rey/11);
84  const scalar turbFrac = 1 - lamFrac;
85 
86  const scalar uStar = sqrt
87  (
88  lamFrac*nuw[facei]*magGradUw[facei] + turbFrac*Cmu5*k[celli]
89  );
90 
91  const scalar omegaVis = 6*nuw[facei]/(beta1_*sqr(y[facei]));
92  const scalar omegaLog = uStar/(Cmu5*nutw.kappa()*y[facei]);
93 
94  G[facei] =
95  lamFrac*Gvf[celli]
96  + turbFrac
97  *sqr(uStar*magGradUw[facei]*y[facei]/uPlus)
98  /(nuw[facei]*nutw.kappa()*yPlus);
99 
100  omega[facei] = lamFrac*omegaVis + turbFrac*omegaLog;
101  }
102  else
103  {
104  if (yPlus < nutw.yPlusLam())
105  {
106  const scalar omegaVis = 6*nuw[facei]/(beta1_*sqr(y[facei]));
107 
108  G[facei] = Gvf[celli];
109 
110  omega[facei] = omegaVis;
111  }
112  else
113  {
114  const scalar uStar = sqrt(Cmu5*k[celli]);
115  const scalar omegaLog = uStar/(Cmu5*nutw.kappa()*y[facei]);
116 
117  G[facei] =
118  sqr(uStar*magGradUw[facei]*y[facei]/uPlus)
119  /(nuw[facei]*nutw.kappa()*yPlus);
120 
121  omega[facei] = omegaLog;
122  }
123  }
124  }
125 
126  return GandOmega;
127 }
128 
129 
130 void Foam::omegaWallFunctionFvPatchScalarField::updateCoeffsMaster()
131 {
132  if (patch().index() != masterPatchIndex())
133  {
134  return;
135  }
136 
137  // Lookup the momentum transport model
138  const momentumTransportModel& mtm =
139  db().lookupType<momentumTransportModel>(internalField().group());
140 
141  // Get mutable references to the turbulence fields
142  VolInternalField<scalar>& G =
143  db().lookupObjectRef<VolInternalField<scalar>>(mtm.GName());
145  const_cast<volScalarField&>
146  (
147  static_cast<const volScalarField&>
148  (
149  internalField()
150  )
151  );
152 
153  const volScalarField::Boundary& omegaBf = omega.boundaryField();
154 
155  // Make all processors build the near wall distances
156  mtm.yb();
157 
158  // Evaluate all the wall functions
159  PtrList<scalarField> Gpfs(omegaBf.size());
160  PtrList<scalarField> omegaPfs(omegaBf.size());
161  forAll(omegaBf, patchi)
162  {
163  if (isA<omegaWallFunctionFvPatchScalarField>(omegaBf[patchi]))
164  {
165  const omegaWallFunctionFvPatchScalarField& omegaPf =
166  refCast<const omegaWallFunctionFvPatchScalarField>
167  (omegaBf[patchi]);
168 
169  Pair<tmp<scalarField>> GandOmega = omegaPf.calculate(mtm);
170 
171  Gpfs.set(patchi, GandOmega.first().ptr());
172  omegaPfs.set(patchi, GandOmega.second().ptr());
173  }
174  }
175 
176  // Average the values into the wall-adjacent cells and store
177  wallCellGPtr_.reset(patchFieldsToWallCellField(Gpfs).ptr());
178  wallCellOmegaPtr_.reset(patchFieldsToWallCellField(omegaPfs).ptr());
179 
180  // Set the fractional proportion of the value in the wall cells
181  UIndirectList<scalar>(G, wallCells()) =
182  (1 - wallCellFraction())*scalarField(G, wallCells())
183  + wallCellFraction()*wallCellGPtr_();
184  UIndirectList<scalar>(omega, wallCells()) =
185  (1 - wallCellFraction())*scalarField(omega, wallCells())
186  + wallCellFraction()*wallCellOmegaPtr_();
187 }
188 
189 
190 void Foam::omegaWallFunctionFvPatchScalarField::manipulateMatrixMaster
191 (
192  fvMatrix<scalar>& matrix
193 )
194 {
195  if (patch().index() != masterPatchIndex())
196  {
197  return;
198  }
199 
200  matrix.setValues(wallCells(), wallCellOmegaPtr_(), wallCellFraction());
201 }
202 
203 
204 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
205 
207 (
208  const fvPatch& p,
210  const dictionary& dict
211 )
212 :
214  beta1_(dict.lookupOrDefault<scalar>("beta1", 0.075)),
215  blended_(dict.lookupOrDefault<Switch>("blended", false)),
216  wallCellGPtr_(nullptr),
217  wallCellOmegaPtr_(nullptr)
218 {}
219 
220 
222 (
224  const fvPatch& p,
226  const fieldMapper& mapper
227 )
228 :
229  wallCellWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
230  beta1_(ptf.beta1_),
231  blended_(ptf.blended_),
232  wallCellGPtr_(nullptr),
233  wallCellOmegaPtr_(nullptr)
234 {}
235 
236 
238 (
241 )
242 :
244  beta1_(owfpsf.beta1_),
245  blended_(owfpsf.blended_),
246  wallCellGPtr_(nullptr),
247  wallCellOmegaPtr_(nullptr)
248 {}
249 
250 
251 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
252 
254 (
255  const fvPatchScalarField& ptf,
256  const fieldMapper& mapper
257 )
258 {
260  wallCellGPtr_.clear();
261  wallCellOmegaPtr_.clear();
262 }
263 
264 
266 (
267  const fvPatchScalarField& ptf
268 )
269 {
271  wallCellGPtr_.clear();
272  wallCellOmegaPtr_.clear();
273 }
274 
275 
277 {
278  if (updated())
279  {
280  return;
281  }
282 
283  initMaster();
284 
285  updateCoeffsMaster();
286 
287  operator==(patchInternalField());
288 
290 }
291 
292 
294 (
295  fvMatrix<scalar>& matrix
296 )
297 {
298  if (manipulatedMatrix())
299  {
300  return;
301  }
302 
303  if (masterPatchIndex() == -1)
304  {
306  << "updateCoeffs must be called before manipulateMatrix"
307  << exit(FatalError);
308  }
309 
310  manipulateMatrixMaster(matrix);
311 
313 }
314 
315 
317 {
318  writeEntry(os, "beta1", beta1_);
319  writeEntry(os, "blended", blended_);
321 }
322 
323 
324 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
325 
326 namespace Foam
327 {
329  (
332  );
333 }
334 
335 
336 // ************************************************************************* //
scalar y
label k
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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...
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
Type of the boundary field.
label size() const
Return the number of elements in the UList.
Definition: ListI.H:171
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:65
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
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
virtual void write(Ostream &) const
Write.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:143
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:202
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
Definition: fvPatchField.C:222
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:356
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:156
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:99
static const nutWallFunctionFvPatchScalarField & nutw(const momentumTransportModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
This boundary condition provides a wall constraint on turbulnce specific dissipation,...
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
omegaWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
Base class for wall functions that modify cell values.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const scalar omega
const scalar yPlus
const scalar uPlus
const scalar Rey
label patchi
compressibleMomentumTransportModel momentumTransportModel
const dimensionedScalar G
Newtonian constant of gravitation.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensionedScalar exp(const dimensionedScalar &ds)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar sqrt(const dimensionedScalar &ds)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
dimensioned< scalar > mag(const dimensioned< Type > &)
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
fvPatchField< vector > fvPatchVectorField
dimensionedScalar pow025(const dimensionedScalar &ds)
dictionary dict
volScalarField & p