omegaWallFunctionFvPatchScalarField.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, 2019 OpenFOAM Foundation
9  Copyright (C) 2017-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 
31 #include "turbulenceModel.H"
32 #include "fvMatrix.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 const Foam::Enum
38 <
39  Foam::omegaWallFunctionFvPatchScalarField::blendingType
40 >
41 Foam::omegaWallFunctionFvPatchScalarField::blendingTypeNames
42 ({
43  { blendingType::STEPWISE , "stepwise" },
44  { blendingType::MAX , "max" },
45  { blendingType::BINOMIAL2 , "binomial2" },
46  { blendingType::BINOMIAL , "binomial" },
47  { blendingType::EXPONENTIAL, "exponential" },
48  { blendingType::TANH, "tanh" }
49 });
50 
52 
53 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
54 
56 {
57  if (master_ != -1)
58  {
59  return;
60  }
61 
62  const volScalarField& omega =
63  static_cast<const volScalarField&>(this->internalField());
64 
65  const volScalarField::Boundary& bf = omega.boundaryField();
66 
67  label master = -1;
68  forAll(bf, patchi)
69  {
70  if (isA<omegaWallFunctionFvPatchScalarField>(bf[patchi]))
71  {
73 
74  if (master == -1)
75  {
76  master = patchi;
77  }
78 
79  opf.master() = master;
80  }
81  }
82 }
83 
84 
86 {
87  const volScalarField& omega =
88  static_cast<const volScalarField&>(this->internalField());
89 
90  const volScalarField::Boundary& bf = omega.boundaryField();
91 
92  const fvMesh& mesh = omega.mesh();
93 
94  if (initialised_ && !mesh.changing())
95  {
96  return;
97  }
98 
99  volScalarField weights
100  (
101  IOobject
102  (
103  "weights",
104  mesh.time().timeName(),
105  mesh,
108  false // do not register
109  ),
110  mesh,
112  );
113 
114  DynamicList<label> omegaPatches(bf.size());
115  forAll(bf, patchi)
116  {
117  if (isA<omegaWallFunctionFvPatchScalarField>(bf[patchi]))
118  {
119  omegaPatches.append(patchi);
120 
121  const labelUList& faceCells = bf[patchi].patch().faceCells();
122  for (const auto& celli : faceCells)
123  {
124  ++weights[celli];
125  }
126  }
127  }
128 
129  cornerWeights_.setSize(bf.size());
130  for (const auto& patchi : omegaPatches)
131  {
132  const fvPatchScalarField& wf = weights.boundaryField()[patchi];
133  cornerWeights_[patchi] = 1.0/wf.patchInternalField();
134  }
135 
136  G_.setSize(internalField().size(), 0.0);
137  omega_.setSize(internalField().size(), 0.0);
138 
139  initialised_ = true;
140 }
141 
142 
145 (
146  const label patchi
147 )
148 {
149  const volScalarField& omega =
150  static_cast<const volScalarField&>(this->internalField());
151 
152  const volScalarField::Boundary& bf = omega.boundaryField();
153 
155  refCast<const omegaWallFunctionFvPatchScalarField>(bf[patchi]);
156 
157  return const_cast<omegaWallFunctionFvPatchScalarField&>(opf);
158 }
159 
160 
162 (
163  const turbulenceModel& turbModel,
164  scalarField& G0,
165  scalarField& omega0
166 )
167 {
168  // accumulate all of the G and omega contributions
169  forAll(cornerWeights_, patchi)
170  {
171  if (!cornerWeights_[patchi].empty())
172  {
173  omegaWallFunctionFvPatchScalarField& opf = omegaPatch(patchi);
174 
175  const List<scalar>& w = cornerWeights_[patchi];
176 
177  opf.calculate(turbModel, w, opf.patch(), G0, omega0);
178  }
179  }
180 
181  // apply zero-gradient condition for omega
182  forAll(cornerWeights_, patchi)
183  {
184  if (!cornerWeights_[patchi].empty())
185  {
186  omegaWallFunctionFvPatchScalarField& opf = omegaPatch(patchi);
187 
188  opf == scalarField(omega0, opf.patch().faceCells());
189  }
190  }
191 }
192 
193 
195 (
196  const turbulenceModel& turbModel,
197  const List<scalar>& cornerWeights,
198  const fvPatch& patch,
199  scalarField& G0,
200  scalarField& omega0
201 )
202 {
203  const label patchi = patch.index();
204 
206  nutWallFunctionFvPatchScalarField::nutw(turbModel, patchi);
207 
208  const scalarField& y = turbModel.y()[patchi];
209 
210  const tmp<scalarField> tnuw = turbModel.nu(patchi);
211  const scalarField& nuw = tnuw();
212 
213  const tmp<volScalarField> tk = turbModel.k();
214  const volScalarField& k = tk();
215 
216  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
217 
218  const scalarField magGradUw(mag(Uw.snGrad()));
219 
220  const scalar Cmu25 = pow025(nutw.Cmu());
221 
222  // Set omega and G
223  forAll(nutw, facei)
224  {
225  const label celli = patch.faceCells()[facei];
226  const scalar yPlus = Cmu25*y[facei]*sqrt(k[celli])/nuw[facei];
227  const scalar w = cornerWeights[facei];
228 
229  // Contribution from the viscous sublayer
230  const scalar omegaVis = 6.0*nuw[facei]/(beta1_*sqr(y[facei]));
231 
232  // Contribution from the inertial sublayer
233  const scalar omegaLog = sqrt(k[celli])/(Cmu25*nutw.kappa()*y[facei]);
234 
235  switch (blending_)
236  {
237  case blendingType::STEPWISE:
238  {
239  if (yPlus > nutw.yPlusLam())
240  {
241  omega0[celli] += w*omegaLog;
242  }
243  else
244  {
245  omega0[celli] += w*omegaVis;
246  }
247  break;
248  }
249 
250  case blendingType::MAX:
251  {
252  // (PH:Eq. 27)
253  omega0[celli] += max(omegaVis, omegaLog);
254  break;
255  }
256 
257  case blendingType::BINOMIAL2:
258  {
259  // (ME:Eq. 15)
260  omega0[celli] += w*sqrt(sqr(omegaVis) + sqr(omegaLog));
261  break;
262  }
263 
264  case blendingType::BINOMIAL:
265  {
266  omega0[celli] +=
267  w*pow
268  (
269  pow(omegaVis, n_) + pow(omegaLog, n_),
270  1.0/n_
271  );
272  break;
273  }
274 
275  case blendingType::EXPONENTIAL:
276  {
277  // (PH:Eq. 31)
278  const scalar Gamma = 0.01*pow4(yPlus)/(1.0 + 5.0*yPlus);
279  const scalar invGamma = 1.0/(Gamma + ROOTVSMALL);
280 
281  omega0[celli] +=
282  w*(omegaVis*exp(-Gamma) + omegaLog*exp(-invGamma));
283  break;
284  }
285 
286  case blendingType::TANH:
287  {
288  // (KAS:Eqs. 33-34)
289  const scalar phiTanh = tanh(pow4(yPlus/10.0));
290  const scalar omegab1 = omegaVis + omegaLog;
291  const scalar omegab2 =
292  pow(pow(omegaVis, 1.2) + pow(omegaLog, 1.2), 1.0/1.2);
293 
294  omega0[celli] += phiTanh*omegab1 + (1.0 - phiTanh)*omegab2;
295  break;
296  }
297  }
298 
299  if (!(blending_ == blendingType::STEPWISE) || yPlus > nutw.yPlusLam())
300  {
301  G0[celli] +=
302  w
303  *(nutw[facei] + nuw[facei])
304  *magGradUw[facei]
305  *Cmu25*sqrt(k[celli])
306  /(nutw.kappa()*y[facei]);
307  }
308  }
309 }
310 
311 
312 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
313 
315 (
316  const fvPatch& p,
318 )
319 :
321  blending_(blendingType::BINOMIAL2),
322  n_(2.0),
323  initialised_(false),
324  master_(-1),
325  beta1_(0.075),
326  G_(),
327  omega_(),
328  cornerWeights_()
329 {}
330 
331 
333 (
335  const fvPatch& p,
337  const fvPatchFieldMapper& mapper
338 )
339 :
340  fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
341  blending_(ptf.blending_),
342  n_(ptf.n_),
343  initialised_(false),
344  master_(-1),
345  beta1_(ptf.beta1_),
346  G_(),
347  omega_(),
348  cornerWeights_()
349 {}
350 
351 
353 (
354  const fvPatch& p,
356  const dictionary& dict
357 )
358 :
360  blending_
361  (
362  blendingTypeNames.getOrDefault
363  (
364  "blending",
365  dict,
366  blendingType::BINOMIAL2
367  )
368  ),
369  n_
370  (
371  dict.getCheckOrDefault<scalar>
372  (
373  "n",
374  2.0,
376  )
377  ),
378  initialised_(false),
379  master_(-1),
380  beta1_(dict.getOrDefault<scalar>("beta1", 0.075)),
381  G_(),
382  omega_(),
383  cornerWeights_()
384 {
385  // The deprecated 'blended' keyword is superseded by the enum 'blending'
386  if (dict.found("blended"))
387  {
389  << "Using deprecated 'blended' keyword"
390  << nl << " Please use either of the below for the same behaviour:"
391  << nl << " 'blending binomial2;' for 'blended on;'"
392  << nl << " 'blending stepwise;' for 'blended off;'"
393  << nl << " OVERWRITING: 'blended' keyword -> 'blending' enum"
394  << endl;
395 
396  bool blended = dict.get<bool>("blended");
397 
398  if (blended)
399  {
400  blending_ = blendingType::BINOMIAL2;
401  }
402  else
403  {
404  blending_ = blendingType::STEPWISE;
405  }
406  }
407 
408  // apply zero-gradient condition on start-up
409  this->operator==(patchInternalField());
410 }
411 
412 
414 (
416 )
417 :
419  blending_(owfpsf.blending_),
420  n_(owfpsf.n_),
421  initialised_(false),
422  master_(-1),
423  beta1_(owfpsf.beta1_),
424  G_(),
425  omega_(),
426  cornerWeights_()
427 {}
428 
429 
431 (
434 )
435 :
436  fixedValueFvPatchField<scalar>(owfpsf, iF),
437  blending_(owfpsf.blending_),
438  n_(owfpsf.n_),
439  initialised_(false),
440  master_(-1),
441  beta1_(owfpsf.beta1_),
442  G_(),
443  omega_(),
444  cornerWeights_()
445 {}
446 
447 
448 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
449 
451 (
452  bool init
453 )
454 {
455  if (patch().index() == master_)
456  {
457  if (init)
458  {
459  G_ = 0.0;
460  }
461 
462  return G_;
463  }
464 
465  return omegaPatch(master_).G();
466 }
467 
468 
470 (
471  bool init
472 )
473 {
474  if (patch().index() == master_)
475  {
476  if (init)
477  {
478  omega_ = 0.0;
479  }
480 
481  return omega_;
482  }
483 
484  return omegaPatch(master_).omega(init);
485 }
486 
487 
489 {
490  if (updated())
491  {
492  return;
493  }
494 
495  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
496  (
498  (
500  internalField().group()
501  )
502  );
503 
504  setMaster();
505 
506  if (patch().index() == master_)
507  {
508  createAveragingWeights();
509  calculateTurbulenceFields(turbModel, G(true), omega(true));
510  }
511 
512  const scalarField& G0 = this->G();
513  const scalarField& omega0 = this->omega();
514 
515  typedef DimensionedField<scalar, volMesh> FieldType;
516 
517  FieldType& G = db().lookupObjectRef<FieldType>(turbModel.GName());
518 
519  FieldType& omega = const_cast<FieldType&>(internalField());
520 
521  forAll(*this, facei)
522  {
523  const label celli = patch().faceCells()[facei];
524 
525  G[celli] = G0[celli];
526  omega[celli] = omega0[celli];
527  }
528 
530 }
531 
532 
534 (
535  const scalarField& weights
536 )
537 {
538  if (updated())
539  {
540  return;
541  }
542 
543  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
544  (
546  (
548  internalField().group()
549  )
550  );
551 
552  setMaster();
553 
554  if (patch().index() == master_)
555  {
556  createAveragingWeights();
557  calculateTurbulenceFields(turbModel, G(true), omega(true));
558  }
559 
560  const scalarField& G0 = this->G();
561  const scalarField& omega0 = this->omega();
562 
563  typedef DimensionedField<scalar, volMesh> FieldType;
564 
565  FieldType& G = db().lookupObjectRef<FieldType>(turbModel.GName());
566 
567  FieldType& omega = const_cast<FieldType&>(internalField());
568 
569  scalarField& omegaf = *this;
570 
571  // only set the values if the weights are > tolerance
572  forAll(weights, facei)
573  {
574  const scalar w = weights[facei];
575 
576  if (w > tolerance_)
577  {
578  const label celli = patch().faceCells()[facei];
579 
580  G[celli] = (1.0 - w)*G[celli] + w*G0[celli];
581  omega[celli] = (1.0 - w)*omega[celli] + w*omega0[celli];
582  omegaf[facei] = omega[celli];
583  }
584  }
585 
587 }
588 
589 
591 (
592  fvMatrix<scalar>& matrix
593 )
594 {
595  if (manipulatedMatrix())
596  {
597  return;
598  }
599 
600  matrix.setValues(patch().faceCells(), patchInternalField());
601 
603 }
604 
605 
607 (
608  fvMatrix<scalar>& matrix,
609  const Field<scalar>& weights
610 )
611 {
612  if (manipulatedMatrix())
613  {
614  return;
615  }
616 
617  DynamicList<label> constraintCells(weights.size());
618  DynamicList<scalar> constraintValues(weights.size());
619  const labelUList& faceCells = patch().faceCells();
620 
621  const DimensionedField<scalar, volMesh>& fld = internalField();
622 
623  forAll(weights, facei)
624  {
625  // only set the values if the weights are > tolerance
626  if (tolerance_ < weights[facei])
627  {
628  const label celli = faceCells[facei];
629 
630  constraintCells.append(celli);
631  constraintValues.append(fld[celli]);
632  }
633  }
634 
635  if (debug)
636  {
637  Pout<< "Patch: " << patch().name()
638  << ": number of constrained cells = " << constraintCells.size()
639  << " out of " << patch().size()
640  << endl;
641  }
642 
643  matrix.setValues(constraintCells, constraintValues);
644 
646 }
647 
648 
650 (
651  Ostream& os
652 ) const
653 {
654  os.writeEntry("blending", blendingTypeNames[blending_]);
655  os.writeEntry("n", n_);
656  os.writeEntry("beta1", beta1_);
658 }
659 
660 
661 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
662 
663 namespace Foam
664 {
666  (
669  );
670 }
671 
672 
673 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::omegaWallFunctionFvPatchScalarField::createAveragingWeights
virtual void createAveragingWeights()
Definition: omegaWallFunctionFvPatchScalarField.C:85
Foam::fvPatchField< scalar >
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField
omegaWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: omegaWallFunctionFvPatchScalarField.C:315
Foam::fvPatchField< scalar >::internalField
const DimensionedField< scalar, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:363
Foam::fvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:225
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::omegaWallFunctionFvPatchScalarField::manipulateMatrix
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
Definition: omegaWallFunctionFvPatchScalarField.C:591
Foam::omegaWallFunctionFvPatchScalarField::calculate
virtual void calculate(const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &omega)
Calculate the omega and G.
Definition: omegaWallFunctionFvPatchScalarField.C:195
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::turbulenceModel::nu
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
Foam::omegaWallFunctionFvPatchScalarField
This boundary condition provides a wall constraint on the specific dissipation rate,...
Definition: omegaWallFunctionFvPatchScalarField.H:291
Foam::nutWallFunctionFvPatchScalarField::Cmu
scalar Cmu() const
Return Cmu.
Definition: nutWallFunctionFvPatchScalarField.H:369
Foam::omegaWallFunctionFvPatchScalarField::updateWeightedCoeffs
virtual void updateWeightedCoeffs(const scalarField &weights)
Update the coefficients associated with the patch field.
Definition: omegaWallFunctionFvPatchScalarField.C:534
Foam::constant::atomic::group
constexpr const char *const group
Group name for atomic constants.
Definition: atomicConstants.H:52
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::DynamicList< label >
Foam::dictionary::found
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Foam::omegaWallFunctionFvPatchScalarField::beta1_
scalar beta1_
beta1 coefficient
Definition: omegaWallFunctionFvPatchScalarField.H:340
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
nutWallFunctionFvPatchScalarField.H
fvMatrix.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:107
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::omegaWallFunctionFvPatchScalarField::master_
label master_
Master patch ID.
Definition: omegaWallFunctionFvPatchScalarField.H:337
Foam::turbulenceModel::y
const nearWallDist & y() const
Return the near wall distances.
Definition: turbulenceModel.H:159
Foam::fixedValueFvPatchField< scalar >
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::fixedValueFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fixedValueFvPatchField.C:156
Foam::pow025
dimensionedScalar pow025(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:133
Foam::omegaWallFunctionFvPatchScalarField::omega
scalarField & omega(bool init=false)
Return non-const access to the master's omega field.
Definition: omegaWallFunctionFvPatchScalarField.C:470
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
Foam::blended
linear/upwind blended differencing scheme.
Definition: blended.H:56
Foam::Field< scalar >
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::tanh
dimensionedScalar tanh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:272
Foam::constant::electromagnetic::G0
const dimensionedScalar G0
Conductance quantum: default SI units: [S].
Foam::omegaWallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: omegaWallFunctionFvPatchScalarField.C:650
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
init
mesh init(true)
Foam::fvPatchField::patchInternalField
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:233
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::omegaWallFunctionFvPatchScalarField::omegaPatch
virtual omegaWallFunctionFvPatchScalarField & omegaPatch(const label patchi)
Helper function to return non-const access to an omega patch.
Definition: omegaWallFunctionFvPatchScalarField.C:145
Foam::dictionary::getCheckOrDefault
T getCheckOrDefault(const word &keyword, const T &deflt, const Predicate &pred, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:209
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
fld
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputLagrangian.H:23
Foam::turbulenceModel::GName
word GName() const
Helper function to return the name of the turbulence G field.
Definition: turbulenceModel.H:138
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:63
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::omegaWallFunctionFvPatchScalarField::calculateTurbulenceFields
virtual void calculateTurbulenceFields(const turbulenceModel &turbulence, scalarField &G0, scalarField &omega0)
Main driver to calculate the turbulence fields.
Definition: omegaWallFunctionFvPatchScalarField.C:162
Foam::turbulenceModel::k
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::nutWallFunctionFvPatchScalarField::yPlusLam
static scalar yPlusLam(const scalar kappa, const scalar E)
Estimate the y+ at the intersection of the two sublayers.
Definition: nutWallFunctionFvPatchScalarField.C:250
Foam::fvPatchField::manipulateMatrix
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
Definition: fvPatchField.C:354
Foam::fvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:321
Foam::fvMatrix::setValues
void setValues(const labelUList &cellLabels, const Type &value)
Definition: fvMatrix.C:1059
Foam::polyMesh::changing
bool changing() const noexcept
Is mesh changing (topology changing and/or moving)
Definition: polyMesh.H:550
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::MinMax::ge
static MinMax< T > ge(const T &minVal)
A semi-infinite range from minVal to the type max.
Definition: MinMaxI.H:31
Foam::nutWallFunctionFvPatchScalarField::nutw
static const nutWallFunctionFvPatchScalarField & nutw(const turbulenceModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
Definition: nutWallFunctionFvPatchScalarField.C:235
Foam::List< scalar >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::nutWallFunctionFvPatchScalarField::kappa
scalar kappa() const
Return kappa.
Definition: nutWallFunctionFvPatchScalarField.H:375
Foam::UList< label >
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::omegaWallFunctionFvPatchScalarField::setMaster
virtual void setMaster()
Definition: omegaWallFunctionFvPatchScalarField.C:55
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
Foam::omegaWallFunctionFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: omegaWallFunctionFvPatchScalarField.C:488
Foam::fvPatch::faceCells
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:113
Foam::omegaWallFunctionFvPatchScalarField::master
virtual label & master()
Return non-const access to the master patch ID.
Definition: omegaWallFunctionFvPatchScalarField.H:387
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
yPlus
scalar yPlus
Definition: evaluateNearWall.H:16
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::fvPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:357
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
IOWarningInFunction
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
Definition: messageStream.H:340
Foam::omegaWallFunctionFvPatchScalarField::tolerance_
static scalar tolerance_
Tolerance used in weighted calculations.
Definition: omegaWallFunctionFvPatchScalarField.H:327
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::nutWallFunctionFvPatchScalarField
The class nutWallFunction is a base class that parents the derived boundary conditions which provide ...
Definition: nutWallFunctionFvPatchScalarField.H:251
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::faceCells
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:56
turbulenceModel.H
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
Foam::turbulenceModel::U
const volVectorField & U() const
Access function to velocity field.
Definition: turbulenceModel.H:144
omegaWallFunctionFvPatchScalarField.H
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::omegaWallFunctionFvPatchScalarField::G
scalarField & G(bool init=false)
Return non-const access to the master's G field.
Definition: omegaWallFunctionFvPatchScalarField.C:451