DAkOmegaSSTFIML.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6  This file is modified from OpenFOAM's source code
7  src/TurbulenceModels/turbulenceModels/RAS/kOmega/kOmega.C
8 
9  OpenFOAM: The Open Source CFD Toolbox
10 
11  Copyright (C): 2011-2016 OpenFOAM Foundation
12 
13  OpenFOAM License:
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 \*---------------------------------------------------------------------------*/
29 
30 #include "DAkOmegaSSTFIML.H"
31 #include "IFstream.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 defineTypeNameAndDebug(DAkOmegaSSTFIML, 0);
39 addToRunTimeSelectionTable(DATurbulenceModel, DAkOmegaSSTFIML, dictionary);
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
43  const word modelType,
44  const fvMesh& mesh,
45  const DAOption& daOption)
46  : DATurbulenceModel(modelType, mesh, daOption),
47  // SST parameters
48  alphaK1_(dimensioned<scalar>::lookupOrAddToDict(
49  "alphaK1",
50  this->coeffDict_,
51  0.85)),
52  alphaK2_(dimensioned<scalar>::lookupOrAddToDict(
53  "alphaK2",
54  this->coeffDict_,
55  1.0)),
56  alphaOmega1_(dimensioned<scalar>::lookupOrAddToDict(
57  "alphaOmega1",
58  this->coeffDict_,
59  0.5)),
60  alphaOmega2_(dimensioned<scalar>::lookupOrAddToDict(
61  "alphaOmega2",
62  this->coeffDict_,
63  0.856)),
64  gamma1_(dimensioned<scalar>::lookupOrAddToDict(
65  "gamma1",
66  this->coeffDict_,
67  5.0 / 9.0)),
68  gamma2_(dimensioned<scalar>::lookupOrAddToDict(
69  "gamma2",
70  this->coeffDict_,
71  0.44)),
72  beta1_(dimensioned<scalar>::lookupOrAddToDict(
73  "beta1",
74  this->coeffDict_,
75  0.075)),
76  beta2_(dimensioned<scalar>::lookupOrAddToDict(
77  "beta2",
78  this->coeffDict_,
79  0.0828)),
80  betaStar_(dimensioned<scalar>::lookupOrAddToDict(
81  "betaStar",
82  this->coeffDict_,
83  0.09)),
84  a1_(dimensioned<scalar>::lookupOrAddToDict(
85  "a1",
86  this->coeffDict_,
87  0.31)),
88  b1_(dimensioned<scalar>::lookupOrAddToDict(
89  "b1",
90  this->coeffDict_,
91  1.0)),
92  c1_(dimensioned<scalar>::lookupOrAddToDict(
93  "c1",
94  this->coeffDict_,
95  10.0)),
96  F3_(Switch::lookupOrAddToDict(
97  "F3",
98  this->coeffDict_,
99  false)),
100  // Augmented variables
101  omega_(const_cast<volScalarField&>(
102  mesh_.thisDb().lookupObject<volScalarField>("omega"))),
103  omegaRes_(
104  IOobject(
105  "omegaRes",
106  mesh.time().timeName(),
107  mesh,
108  IOobject::NO_READ,
109  IOobject::NO_WRITE),
110  mesh,
111 #ifdef CompressibleFlow
112  dimensionedScalar("omegaRes", dimensionSet(1, -3, -2, 0, 0, 0, 0), 0.0),
113 #endif
114 #ifdef IncompressibleFlow
115  dimensionedScalar("omegaRes", dimensionSet(0, 0, -2, 0, 0, 0, 0), 0.0),
116 #endif
117  zeroGradientFvPatchField<scalar>::typeName),
118  omegaResRef_(
119  IOobject(
120  "omegaResRef",
121  mesh.time().timeName(),
122  mesh,
123  IOobject::NO_READ,
124  IOobject::NO_WRITE),
125  omegaRes_),
126  omegaResPartDeriv_(
127  IOobject(
128  "omegaResPartDeriv",
129  mesh.time().timeName(),
130  mesh,
131  IOobject::NO_READ,
132  IOobject::NO_WRITE),
133  omegaRes_),
134  omegaRef_(
135  IOobject(
136  "omegaRef",
137  mesh.time().timeName(),
138  mesh,
139  IOobject::NO_READ,
140  IOobject::NO_WRITE),
141  omega_),
142  k_(
143  const_cast<volScalarField&>(
144  mesh_.thisDb().lookupObject<volScalarField>("k"))),
145  kRes_(
146  IOobject(
147  "kRes",
148  mesh.time().timeName(),
149  mesh,
150  IOobject::NO_READ,
151  IOobject::NO_WRITE),
152  mesh,
153 #ifdef CompressibleFlow
154  dimensionedScalar("kRes", dimensionSet(1, -1, -3, 0, 0, 0, 0), 0.0),
155 #endif
156 #ifdef IncompressibleFlow
157  dimensionedScalar("kRes", dimensionSet(0, 2, -3, 0, 0, 0, 0), 0.0),
158 #endif
159  zeroGradientFvPatchField<scalar>::typeName),
160  kResRef_(
161  IOobject(
162  "kResRef",
163  mesh.time().timeName(),
164  mesh,
165  IOobject::NO_READ,
166  IOobject::NO_WRITE),
167  kRes_),
168  kResPartDeriv_(
169  IOobject(
170  "kResPartDeriv",
171  mesh.time().timeName(),
172  mesh,
173  IOobject::NO_READ,
174  IOobject::NO_WRITE),
175  kRes_),
176  kRef_(
177  IOobject(
178  "kRef",
179  mesh.time().timeName(),
180  mesh,
181  IOobject::NO_READ,
182  IOobject::NO_WRITE),
183  k_),
185  betaFieldInversion_(
186  IOobject(
187  "betaFieldInversion",
188  mesh.time().timeName(),
189  mesh_,
190  IOobject::READ_IF_PRESENT,
191  IOobject::AUTO_WRITE),
192  mesh_,
193  dimensionedScalar("betaFieldInversion", dimensionSet(0, 0, 0, 0, 0, 0, 0), 1.0),
194  zeroGradientFvPatchScalarField::typeName),
195  betaFieldInversionML_(
196  IOobject(
197  "betaFieldInversionML",
198  mesh.time().timeName(),
199  mesh_,
200  IOobject::NO_READ,
201  IOobject::AUTO_WRITE),
202  mesh_,
203  dimensionedScalar("betaFieldInversionML", dimensionSet(0, 0, 0, 0, 0, 0, 0), 1.0),
204  zeroGradientFvPatchScalarField::typeName),
205  QCriterion_(
206  IOobject(
207  "QCriterion",
208  mesh.time().timeName(),
209  mesh_,
210  IOobject::NO_READ,
211  IOobject::AUTO_WRITE),
212  mesh_,
213  dimensionedScalar("QCriterion", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0),
214  zeroGradientFvPatchScalarField::typeName),
215  p_(const_cast<volScalarField&>(
216  mesh_.thisDb().lookupObject<volScalarField>("p"))),
217  pGradAlongStream_(
218  IOobject(
219  "pGradAlongStream",
220  mesh.time().timeName(),
221  mesh_,
222  IOobject::NO_READ,
223  IOobject::AUTO_WRITE),
224  mesh_,
225  dimensionedScalar("pressureAlongStream", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0),
226  zeroGradientFvPatchScalarField::typeName),
227  turbulenceIntensity_(
228  IOobject(
229  "turbulenceIntensity",
230  mesh.time().timeName(),
231  mesh_,
232  IOobject::NO_READ,
233  IOobject::AUTO_WRITE),
234  mesh_,
235  dimensionedScalar("turbulenceIntensity", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0),
236  zeroGradientFvPatchScalarField::typeName),
237  transportProperties_(
238  IOobject(
239  "transportProperties",
240  mesh.time().constant(),
241  mesh_,
242  IOobject::MUST_READ,
243  IOobject::NO_WRITE)),
244  ReT_(
245  IOobject(
246  "ReT",
247  mesh.time().timeName(),
248  mesh_,
249  IOobject::NO_READ,
250  IOobject::AUTO_WRITE),
251  mesh_,
252  dimensionedScalar("ReT", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0),
253  zeroGradientFvPatchScalarField::typeName),
254  convectionTKE_(
255  IOobject(
256  "convectionTKE",
257  mesh.time().timeName(),
258  mesh_,
259  IOobject::NO_READ,
260  IOobject::AUTO_WRITE),
261  mesh_,
262  dimensionedScalar("convectionTKE", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0),
263  zeroGradientFvPatchScalarField::typeName),
264  tauRatio_(
265  IOobject(
266  "tauRatio",
267  mesh.time().timeName(),
268  mesh_,
269  IOobject::NO_READ,
270  IOobject::AUTO_WRITE),
271  mesh_,
272  dimensionedScalar("tauRatio", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0),
273  zeroGradientFvPatchScalarField::typeName),
274  pressureStress_(
275  IOobject(
276  "pressureStress",
277  mesh.time().timeName(),
278  mesh_,
279  IOobject::NO_READ,
280  IOobject::AUTO_WRITE),
281  mesh_,
282  dimensionedScalar("pressureStress", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0),
283  zeroGradientFvPatchScalarField::typeName),
284  curvature_(
285  IOobject(
286  "curvature",
287  mesh.time().timeName(),
288  mesh_,
289  IOobject::NO_READ,
290  IOobject::AUTO_WRITE),
291  mesh_,
292  dimensionedScalar("curvature", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0),
293  zeroGradientFvPatchScalarField::typeName),
294  UGradMisalignment_(
295  IOobject(
296  "UGradMisalignment",
297  mesh.time().timeName(),
298  mesh_,
299  IOobject::NO_READ,
300  IOobject::AUTO_WRITE),
301  mesh_,
302  dimensionedScalar("UGradMisalignment", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0),
303  zeroGradientFvPatchScalarField::typeName),
304  y_(mesh_.thisDb().lookupObject<volScalarField>("yWall"))
305 {
306 
307  // initialize printInterval_ we need to check whether it is a steady state
308  // or unsteady primal solver
309  IOdictionary fvSchemes(
310  IOobject(
311  "fvSchemes",
312  mesh.time().system(),
313  mesh,
314  IOobject::MUST_READ,
315  IOobject::NO_WRITE,
316  false));
317  word ddtScheme = word(fvSchemes.subDict("ddtSchemes").lookup("default"));
318  if (ddtScheme == "steadyState")
319  {
321  daOption.getAllOptions().lookupOrDefault<label>("printInterval", 100);
322  }
323  else
324  {
326  daOption.getAllOptions().lookupOrDefault<label>("printIntervalUnsteady", 500);
327  }
328 
329  // calculate the size of omegaWallFunction faces
330  label nWallFaces = 0;
331  forAll(omega_.boundaryField(), patchI)
332  {
333  if (omega_.boundaryField()[patchI].type() == "omegaWallFunction"
334  && omega_.boundaryField()[patchI].size() > 0)
335  {
336  forAll(omega_.boundaryField()[patchI], faceI)
337  {
338  nWallFaces++;
339  }
340  }
341  }
342 
343  // initialize omegaNearWall
344  omegaNearWall_.setSize(nWallFaces);
345 }
346 
347 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
348 
349 // SA member functions. these functions are copied from
350 tmp<volScalarField> DAkOmegaSSTFIML::F1(
351  const volScalarField& CDkOmega) const
352 {
353 
354  tmp<volScalarField> CDkOmegaPlus = max(
355  CDkOmega,
356  dimensionedScalar("1.0e-10", dimless / sqr(dimTime), 1.0e-10));
357 
358  tmp<volScalarField> arg1 = min(
359  min(
360  max(
361  (scalar(1) / betaStar_) * sqrt(k_) / (omega_ * y_),
362  scalar(500) * (this->nu()) / (sqr(y_) * omega_)),
363  (scalar(4) * alphaOmega2_) * k_ / (CDkOmegaPlus * sqr(y_))),
364  scalar(10));
365 
366  return tanh(pow4(arg1));
367 }
368 
369 tmp<volScalarField> DAkOmegaSSTFIML::F2() const
370 {
371 
372  tmp<volScalarField> arg2 = min(
373  max(
374  (scalar(2) / betaStar_) * sqrt(k_) / (omega_ * y_),
375  scalar(500) * (this->nu()) / (sqr(y_) * omega_)),
376  scalar(100));
377 
378  return tanh(sqr(arg2));
379 }
380 
381 tmp<volScalarField> DAkOmegaSSTFIML::F3() const
382 {
383 
384  tmp<volScalarField> arg3 = min(
385  150 * (this->nu()) / (omega_ * sqr(y_)),
386  scalar(10));
387 
388  return 1 - tanh(pow4(arg3));
389 }
390 
391 tmp<volScalarField> DAkOmegaSSTFIML::F23() const
392 {
393  tmp<volScalarField> f23(F2());
394 
395  if (F3_)
396  {
397  f23.ref() *= F3();
398  }
399 
400  return f23;
401 }
402 
403 tmp<volScalarField::Internal> DAkOmegaSSTFIML::GbyNu(
404  const volScalarField::Internal& GbyNu0,
405  const volScalarField::Internal& F2,
406  const volScalarField::Internal& S2) const
407 {
408  return min(
409  GbyNu0,
410  (c1_ / a1_) * betaStar_ * omega_()
411  * max(a1_ * omega_(), b1_ * F2 * sqrt(S2)));
412 }
413 
414 tmp<volScalarField::Internal> DAkOmegaSSTFIML::Pk(
415  const volScalarField::Internal& G) const
416 {
417  return min(G, (c1_ * betaStar_) * k_() * omega_());
418 }
419 
420 tmp<volScalarField::Internal> DAkOmegaSSTFIML::epsilonByk(
421  const volScalarField& F1,
422  const volTensorField& gradU) const
423 {
424  return betaStar_ * omega_();
425 }
426 
427 tmp<fvScalarMatrix> DAkOmegaSSTFIML::kSource() const
428 {
429  const volScalarField& rho = rho_;
430  return tmp<fvScalarMatrix>(
431  new fvScalarMatrix(
432  k_,
433  dimVolume * rho.dimensions() * k_.dimensions() / dimTime));
434 }
435 
436 tmp<fvScalarMatrix> DAkOmegaSSTFIML::omegaSource() const
437 {
438  const volScalarField& rho = rho_;
439  return tmp<fvScalarMatrix>(
440  new fvScalarMatrix(
441  omega_,
442  dimVolume * rho.dimensions() * omega_.dimensions() / dimTime));
443 }
444 
445 tmp<fvScalarMatrix> DAkOmegaSSTFIML::Qsas(
446  const volScalarField::Internal& S2,
447  const volScalarField::Internal& gamma,
448  const volScalarField::Internal& beta) const
449 {
450  const volScalarField& rho = rho_;
451  return tmp<fvScalarMatrix>(
452  new fvScalarMatrix(
453  omega_,
454  dimVolume * rho.dimensions() * omega_.dimensions() / dimTime));
455 }
456 
457 // Augmented functions
458 void DAkOmegaSSTFIML::correctModelStates(wordList& modelStates) const
459 {
460  /*
461  Description:
462  Update the name in modelStates based on the selected physical model at runtime
463 
464  Example:
465  In DAStateInfo, if the modelStates reads:
466 
467  modelStates = {"nut"}
468 
469  then for the SA model, calling correctModelStates(modelStates) will give:
470 
471  modelStates={"nuTilda"}
472 
473  while calling correctModelStates(modelStates) for the SST model will give
474 
475  modelStates={"k","omega"}
476 
477  We don't udpate the names for the radiation model because users are
478  supposed to set modelStates={"G"}
479  */
480 
481  // For SST model, we need to replace nut with k, omega
482 
483  forAll(modelStates, idxI)
484  {
485  word stateName = modelStates[idxI];
486  if (stateName == "nut")
487  {
488  modelStates[idxI] = "omega";
489  modelStates.append("k");
490  }
491  }
492 }
493 
495 {
496  /*
497  Description:
498  Update nut based on other turbulence variables and update the BCs
499  Also update alphat if is present
500  */
501 
502  const volVectorField U = mesh_.thisDb().lookupObject<volVectorField>("U");
503  tmp<volTensorField> tgradU = fvc::grad(U);
504  volScalarField S2(2 * magSqr(symm(tgradU())));
505 
506  nut_ = a1_ * k_ / max(a1_ * omega_, b1_ * F23() * sqrt(S2));
507 
508  nut_.correctBoundaryConditions(); // nutkWallFunction: update wall face nut based on k
509 
510  // this is basically BasicTurbulenceModel::correctNut();
511  this->correctAlphat();
512 
513  return;
514 }
515 
517 {
518  /*
519  Description:
520  Update turbulence variable boundary values
521  */
522 
523  // correct the BCs for the perturbed fields
524  // kqWallFunction is a zero-gradient BC
525  k_.correctBoundaryConditions();
526 }
527 
529 {
530  /*
531  Description:
532  this is a special treatment for omega BC because we cant directly call omega.
533  correctBoundaryConditions() because it will modify the internal omega and G that
534  are right next to walls. This will mess up adjoint Jacobians
535  To solve this issue,
536  1. we store the near wall omega before calling omega.correctBoundaryConditions()
537  2. call omega.correctBoundaryConditions()
538  3. Assign the stored near wall omega back
539  4. Apply a zeroGradient BC for omega at the wall patches
540  *********** NOTE *************
541  this treatment will obviously downgrade the accuracy of adjoint derivative since it is
542  not 100% consistent with what is used for the flow solver; however, our observation
543  shows that the impact is not very large.
544  */
545 
546  // save the perturbed omega at the wall
547  this->saveOmegaNearWall();
548  // correct omega boundary conditions, this includes updating wall face and near wall omega values,
549  // updating the inter-proc BCs
550  omega_.correctBoundaryConditions();
551  // reset the corrected omega near wall cell to its perturbed value
552  this->setOmegaNearWall();
553 }
554 
556 {
557  /*
558  Description:
559  Save the current near wall omega values to omegaNearWall_
560  */
561  label counterI = 0;
562  forAll(omega_.boundaryField(), patchI)
563  {
564  if (omega_.boundaryField()[patchI].type() == "omegaWallFunction"
565  and omega_.boundaryField()[patchI].size() > 0)
566  {
567  const UList<label>& faceCells = mesh_.boundaryMesh()[patchI].faceCells();
568  forAll(faceCells, faceI)
569  {
570  //Info<<"Near Wall cellI: "<<faceCells[faceI]<<endl;
571  omegaNearWall_[counterI] = omega_[faceCells[faceI]];
572  counterI++;
573  }
574  }
575  }
576  return;
577 }
578 
580 {
581  /*
582  Description:
583  Set the current near wall omega values from omegaNearWall_
584  Here we also apply a zeroGradient BC to the wall faces
585  */
586  label counterI = 0;
587  forAll(omega_.boundaryField(), patchI)
588  {
589  if (omega_.boundaryField()[patchI].type() == "omegaWallFunction"
590  && omega_.boundaryField()[patchI].size() > 0)
591  {
592  const UList<label>& faceCells = mesh_.boundaryMesh()[patchI].faceCells();
593  forAll(faceCells, faceI)
594  {
595  omega_[faceCells[faceI]] = omegaNearWall_[counterI];
596  // zeroGradient BC
597  omega_.boundaryFieldRef()[patchI][faceI] = omega_[faceCells[faceI]];
598  counterI++;
599  }
600  }
601  }
602  return;
603 }
604 
606 {
607  /*
608  Description:
609  Update nut based on nuTilda. Note: we need to update nut and its BC since we
610  may have perturbed other turbulence vars that affect the nut values
611  */
612 
613  this->correctNut();
614 }
615 
616 void DAkOmegaSSTFIML::correctStateResidualModelCon(List<List<word>>& stateCon) const
617 {
618  /*
619  Description:
620  Update the original variable connectivity for the adjoint state
621  residuals in stateCon. Basically, we modify/add state variables based on the
622  original model variables defined in stateCon.
623 
624  Input:
625 
626  stateResCon: the connectivity levels for a state residual, defined in Foam::DAJacCon
627 
628  Example:
629  If stateCon reads:
630  stateCon=
631  {
632  {"U", "p", "nut"},
633  {"p"}
634  }
635 
636  For the SA turbulence model, calling this function for will get a new stateCon
637  stateCon=
638  {
639  {"U", "p", "nuTilda"},
640  {"p"}
641  }
642 
643  For the SST turbulence model, calling this function will give
644  stateCon=
645  {
646  {"U", "p", "k", "omega"},
647  {"p", "U"}
648  }
649  ***NOTE***: we add a extra level of U connectivity because nut is
650  related to grad(U), k, and omega in SST!
651  */
652 
653  label stateConSize = stateCon.size();
654  forAll(stateCon, idxI)
655  {
656  label addUCon = 0;
657  forAll(stateCon[idxI], idxJ)
658  {
659  word conStateName = stateCon[idxI][idxJ];
660  if (conStateName == "nut")
661  {
662  stateCon[idxI][idxJ] = "omega"; // replace nut with omega
663  stateCon[idxI].append("k"); // also add k for that level
664  addUCon = 1;
665  }
666  }
667  // add U for the current level and level+1 if it is not there yet
668  label isU;
669  if (addUCon == 1)
670  {
671  // first add U for the current level
672  isU = 0;
673  forAll(stateCon[idxI], idxJ)
674  {
675  word conStateName = stateCon[idxI][idxJ];
676  if (conStateName == "U")
677  {
678  isU = 1;
679  }
680  }
681  if (!isU)
682  {
683  stateCon[idxI].append("U");
684  }
685 
686  // now add U for level+1 if idxI is not the largest level
687  // if idxI is already the largest level, we have a problem
688  if (idxI != stateConSize - 1)
689  {
690  isU = 0;
691  forAll(stateCon[idxI + 1], idxJ)
692  {
693  word conStateName = stateCon[idxI + 1][idxJ];
694  if (conStateName == "U")
695  {
696  isU = 1;
697  }
698  }
699  if (!isU)
700  {
701  stateCon[idxI + 1].append("U");
702  }
703  }
704  else
705  {
706  FatalErrorIn(
707  "In DAStateInfo, nut shows in the largest connectivity level! "
708  "This is not supported!")
709  << exit(FatalError);
710  }
711  }
712  }
713 }
714 
715 void DAkOmegaSSTFIML::addModelResidualCon(HashTable<List<List<word>>>& allCon) const
716 {
717  /*
718  Description:
719  Add the connectivity levels for all physical model residuals to allCon
720 
721  Input:
722  allCon: the connectivity levels for all state residual, defined in DAJacCon
723 
724  Example:
725  If stateCon reads:
726  allCon=
727  {
728  "URes":
729  {
730  {"U", "p", "nut"},
731  {"p"}
732  }
733  }
734 
735  For the SA turbulence model, calling this function for will get a new stateCon,
736  something like this:
737  allCon=
738  {
739  "URes":
740  {
741  {"U", "p", "nuTilda"},
742  {"p"}
743  },
744  "nuTildaRes":
745  {
746  {"U", "phi", "nuTilda"},
747  {"U"}
748  }
749  }
750 
751  */
752 
753  word pName;
754 
755  if (mesh_.thisDb().foundObject<volScalarField>("p"))
756  {
757  pName = "p";
758  }
759  else if (mesh_.thisDb().foundObject<volScalarField>("p_rgh"))
760  {
761  pName = "p_rgh";
762  }
763  else
764  {
765  FatalErrorIn(
766  "Neither p nor p_rgh was found in mesh.thisDb()!"
767  "addModelResidualCon failed to setup turbulence residuals!")
768  << exit(FatalError);
769  }
770 
771  // NOTE: for compressible flow, it depends on rho so we need to add T and p
772 #ifdef IncompressibleFlow
773  allCon.set(
774  "omegaRes",
775  {
776  {"U", "omega", "k", "phi"}, // lv0
777  {"U", "omega", "k"}, // lv1
778  {"U", "omega", "k"} // lv2
779  });
780  allCon.set(
781  "kRes",
782  {
783  {"U", "omega", "k", "phi"}, // lv0
784  {"U", "omega", "k"}, // lv1
785  {"U", "omega", "k"} // lv2
786  });
787 #endif
788 
789 #ifdef CompressibleFlow
790  allCon.set(
791  "omegaRes",
792  {
793  {"U", "T", pName, "omega", "k", "phi"}, // lv0
794  {"U", "T", pName, "omega", "k"}, // lv1
795  {"U", "T", pName, "omega", "k"} // lv2
796  });
797  allCon.set(
798  "kRes",
799  {
800  {"U", "T", pName, "omega", "k", "phi"}, // lv0
801  {"U", "T", pName, "omega", "k"}, // lv1
802  {"U", "T", pName, "omega", "k"} // lv2
803  });
804 #endif
805 }
806 
808 {
809  /*
810  Descroption:
811  Solve the residual equations and update the state. This function will be called
812  by the DASolver. It is needed because we want to control the output frequency
813  of the residual convergence every 100 steps. If using the correct from turbulence
814  it will output residual every step which will be too much of information.
815  */
816 
817  // We set the flag solveTurbState_ to 1 such that in the calcResiduals function
818  // we will solve and update nuTilda
819  solveTurbState_ = 1;
820  dictionary dummyOptions;
821  this->calcResiduals(dummyOptions);
822  // after it, we reset solveTurbState_ = 0 such that calcResiduals will not
823  // update nuTilda when calling from the adjoint class, i.e., solveAdjoint from DASolver.
824  solveTurbState_ = 0;
825 }
826 
828 {
829 
830  // COMPUTE MACHINE LEARNING FEATURES
832  volTensorField UGrad(fvc::grad(U_));
833  volTensorField Omega("Omega", skew(UGrad));
834  volScalarField magOmegaSqr(magSqr(Omega));
835  volSymmTensorField S("S", symm(UGrad));
836  volScalarField magS(mag(S));
837  volScalarField magSSqr(magSqr(S));
838  QCriterion_ = (magOmegaSqr - magSSqr) / (magOmegaSqr + magSSqr);
839 
841  volVectorField pGrad("gradP", fvc::grad(p_));
842  volScalarField pG_denominator(mag(U_) * mag(pGrad) + mag(U_ & pGrad));
843  pGradAlongStream_ = (U_ & pGrad) / Foam::max(pG_denominator, dimensionedScalar("minpG", dimensionSet(0, 2, -3, 0, 0, 0, 0), SMALL));
844 
846  turbulenceIntensity_ = k_ / (0.5 * (U_ & U_) + k_);
847 
849  //dimensionedScalar refViscosity(this->nu());
850  dimensionedScalar maxReT("maxReT", dimless, 2.0);
851  ReT_ = Foam::min((sqrt(k_) * y_) / (scalar(50.0) * this->nu()), maxReT);
852 
854  volSymmTensorField tau(2.0 / 3.0 * I * k_ - nut_ * twoSymm(fvc::grad(U_)));
855  volVectorField kGrad("gradK", fvc::grad(k_));
856  convectionTKE_ = (U_ & kGrad) / (mag(tau && S) + mag(U_ & kGrad));
857 
859  tauRatio_ = mag(tau) / (k_ + mag(tau));
860 
862  volVectorField diagUGrad(
863  IOobject("diagUGrad",
864  mesh_.time().timeName(),
865  mesh_,
866  IOobject::NO_READ,
867  IOobject::NO_WRITE),
868  mesh_,
869  dimensionedVector("diagUGrad", dimensionSet(0, 0, 0, 0, 0, 0, 0), Foam::vector(0, 0, 0)),
870  zeroGradientFvPatchScalarField::typeName);
871 
872  forAll(mesh_.C(), cI)
873  {
874  diagUGrad[cI].component(0) = UGrad[cI].xx();
875  diagUGrad[cI].component(1) = UGrad[cI].yy();
876  diagUGrad[cI].component(2) = UGrad[cI].zz();
877  pressureStress_[cI] = mag(pGrad[cI]) / (mag(pGrad[cI]) + mag(3.0 * cmptAv(U_[cI] & diagUGrad[cI])));
878  }
879 
881  forAll(mesh_.C(), cI)
882  {
883  curvature_[cI] = mag(U_[cI] & UGrad[cI]) / (mag(U_[cI] & U_[cI]) + mag(U_[cI] & UGrad[cI]));
884  }
885 
887  forAll(mesh_.C(), cI)
888  {
889  UGradMisalignment_[cI] = mag(U_[cI] & UGrad[cI] & U_[cI])
890  / (mag(U_[cI]) * mag(UGrad[cI] & U_[cI]) + mag(U_[cI] & UGrad[cI] & U_[cI]));
891  }
892 
893  label n = 9 * mesh_.nCells();
894  label m = mesh_.nCells();
895 
896  forAll(mesh_.cells(), cI)
897  {
898  inputs_[cI * 9 + 0] = QCriterion_[cI];
899  inputs_[cI * 9 + 1] = UGradMisalignment_[cI];
900  inputs_[cI * 9 + 2] = pGradAlongStream_[cI];
901  inputs_[cI * 9 + 3] = turbulenceIntensity_[cI];
902  inputs_[cI * 9 + 4] = ReT_[cI];
903  inputs_[cI * 9 + 5] = convectionTKE_[cI];
904  inputs_[cI * 9 + 6] = curvature_[cI];
905  inputs_[cI * 9 + 7] = pressureStress_[cI];
906  inputs_[cI * 9 + 8] = tauRatio_[cI];
907  }
908 
909  // NOTE: forward mode not supported..
910 #if defined(CODI_AD_REVERSE)
911 
912  // we need to use the external function helper from CoDiPack to propagate the AD
913 
914  codi::ExternalFunctionHelper<codi::RealReverse> externalFunc;
915  for (label i = 0; i < mesh_.nCells() * 9; i++)
916  {
917  externalFunc.addInput(inputs_[i]);
918  }
919 
920  for (label i = 0; i < mesh_.nCells(); i++)
921  {
922  externalFunc.addOutput(outputs_[i]);
923  }
924 
925  externalFunc.callPrimalFunc(DAkOmegaSSTFIML::betaCompute);
926 
927  codi::RealReverse::Tape& tape = codi::RealReverse::getTape();
928 
929  if (tape.isActive())
930  {
931  externalFunc.addToTape(DAkOmegaSSTFIML::betaJacVecProd);
932  }
933 
935  {
936  betaFieldInversionML_[cellI] = outputs_[cellI];
937  }
938 
939 #elif defined(CODI_AD_FORWARD)
940 
941  for (label i = 0; i < n; i++)
942  {
943  inputsDouble_[i] = inputs_[i].value();
944  }
945 
946  for (label i = 0; i < m; i++)
947  {
948  outputsDouble_[i] = outputs_[i].value();
949  }
950 
951  // python callback function
952  DAUtility::pyCalcBetaInterface(inputsDouble_, n, outputsDouble_, m, DAUtility::pyCalcBeta);
953 
955  {
956  betaFieldInversionML_[cellI] = outputsDouble_[cellI];
957  }
958 
959 #else
960 
961  // python callback function
963 
965  {
966  betaFieldInversionML_[cellI] = outputs_[cellI];
967  }
968 
969 #endif
970 }
971 
972 void DAkOmegaSSTFIML::calcResiduals(const dictionary& options)
973 {
974  /*
975  Descroption:
976  If solveTurbState_ == 1, this function solve and update k and omega, and
977  is the same as calling turbulence.correct(). If solveTurbState_ == 0,
978  this function compute residuals for turbulence variables, e.g., nuTildaRes_
979 
980  Input:
981  options.isPC: 1 means computing residuals for preconditioner matrix.
982  This essentially use the first order scheme for div(phi,nuTilda)
983 
984  p_, U_, phi_, etc: State variables in OpenFOAM
985 
986  Output:
987  kRes_/omegaRes_: If solveTurbState_ == 0, update the residual field variable
988 
989  k_/omega_: If solveTurbState_ == 1, update them
990  */
991 
992  // Copy and modify based on the "correct" function
993 
994  label printToScreen = this->isPrintTime(mesh_.time(), printInterval_);
995 
996  word divKScheme = "div(phi,k)";
997  word divOmegaScheme = "div(phi,omega)";
998 
999  label isPC = 0;
1000 
1001  if (!solveTurbState_)
1002  {
1003  isPC = options.getLabel("isPC");
1004 
1005  if (isPC)
1006  {
1007  divKScheme = "div(pc)";
1008  divOmegaScheme = "div(pc)";
1009  }
1010  }
1011 
1012  this->calcBetaField();
1013 
1014  // Note: for compressible flow, the "this->phi()" function divides phi by fvc:interpolate(rho),
1015  // while for the incompresssible "this->phi()" returns phi only
1016  // see src/TurbulenceModels/compressible/compressibleTurbulenceModel.C line 62 to 73
1017  volScalarField::Internal divU(fvc::div(fvc::absolute(phi_ / fvc::interpolate(rho_), U_)));
1018 
1019  tmp<volTensorField> tgradU = fvc::grad(U_);
1020  volScalarField S2(2 * magSqr(symm(tgradU())));
1021  volScalarField::Internal GbyNu0((tgradU() && dev(twoSymm(tgradU()))));
1022  volScalarField::Internal G("kOmegaSSTFIML:G", nut_ * GbyNu0);
1023 
1024  if (solveTurbState_)
1025  {
1026  // Update omega and G at the wall
1027  omega_.boundaryFieldRef().updateCoeffs();
1028  }
1029  else
1030  {
1031  // NOTE instead of calling omega_.boundaryFieldRef().updateCoeffs();
1032  // here we call our self-defined boundary conditions
1034  }
1035 
1036  volScalarField CDkOmega(
1037  (scalar(2) * alphaOmega2_) * (fvc::grad(k_) & fvc::grad(omega_)) / omega_);
1038 
1039  volScalarField F1(this->F1(CDkOmega));
1040  volScalarField F23(this->F23());
1041 
1042  {
1043 
1044  volScalarField::Internal gamma(this->gamma(F1));
1045  volScalarField::Internal beta(this->beta(F1));
1046 
1047  // Turbulent frequency equation
1048  tmp<fvScalarMatrix> omegaEqn(
1049  fvm::ddt(phase_, rho_, omega_)
1050  + fvm::div(phaseRhoPhi_, omega_, divOmegaScheme)
1051  - fvm::laplacian(phase_ * rho_ * DomegaEff(F1), omega_)
1052  == betaFieldInversionML_() * phase_() * rho_() * gamma * GbyNu(GbyNu0, F23(), S2())
1053  - fvm::SuSp((2.0 / 3.0) * phase_() * rho_() * gamma * divU, omega_)
1054  - fvm::Sp(phase_() * rho_() * beta * omega_(), omega_)
1055  - fvm::SuSp(
1056  phase_() * rho_() * (F1() - scalar(1)) * CDkOmega() / omega_(),
1057  omega_)
1058  + Qsas(S2(), gamma, beta)
1059  + omegaSource()
1060 
1061  );
1062 
1063  omegaEqn.ref().relax();
1064  omegaEqn.ref().boundaryManipulate(omega_.boundaryFieldRef());
1065 
1066  if (solveTurbState_)
1067  {
1068 
1069  // get the solver performance info such as initial
1070  // and final residuals
1071  SolverPerformance<scalar> solverOmega = solve(omegaEqn);
1072  if (printToScreen)
1073  {
1074  Info << "omega Initial residual: " << solverOmega.initialResidual() << endl
1075  << " Final residual: " << solverOmega.finalResidual() << endl;
1076  }
1077 
1078  DAUtility::boundVar(allOptions_, omega_, printToScreen);
1079  }
1080  else
1081  {
1082  // reset the corrected omega near wall cell to its perturbed value
1083  this->setOmegaNearWall();
1084 
1085  // calculate residuals
1086  omegaRes_ = omegaEqn() & omega_;
1087  // need to normalize residuals
1088  normalizeResiduals(omegaRes);
1089  }
1090  }
1091 
1092  // Turbulent kinetic energy equation
1093  tmp<fvScalarMatrix> kEqn(
1094  fvm::ddt(phase_, rho_, k_)
1095  + fvm::div(phaseRhoPhi_, k_, divKScheme)
1096  - fvm::laplacian(phase_ * rho_ * DkEff(F1), k_)
1097  == phase_() * rho_() * Pk(G)
1098  - fvm::SuSp((2.0 / 3.0) * phase_() * rho_() * divU, k_)
1099  - fvm::Sp(phase_() * rho_() * epsilonByk(F1, tgradU()), k_)
1100  + kSource());
1101 
1102  tgradU.clear();
1103 
1104  kEqn.ref().relax();
1105 
1106  if (solveTurbState_)
1107  {
1108 
1109  // get the solver performance info such as initial
1110  // and final residuals
1111  SolverPerformance<scalar> solverK = solve(kEqn);
1112  if (printToScreen)
1113  {
1114  Info << "k Initial residual: " << solverK.initialResidual() << endl
1115  << " Final residual: " << solverK.finalResidual() << endl;
1116  }
1117 
1118  DAUtility::boundVar(allOptions_, k_, printToScreen);
1119 
1120  this->correctNut();
1121  }
1122  else
1123  {
1124  // calculate residuals
1125  kRes_ = kEqn() & k_;
1126  // need to normalize residuals
1127  normalizeResiduals(kRes);
1128  }
1129 
1130  return;
1131 }
1132 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1133 
1134 } // End namespace Foam
1135 
1136 // ************************************************************************* //
Foam::DAkOmegaSSTFIML::DAkOmegaSSTFIML
DAkOmegaSSTFIML(const word modelType, const fvMesh &mesh, const DAOption &daOption)
Definition: DAkOmegaSSTFIML.C:42
Foam::DAkOmegaSSTFIML::F3_
Switch F3_
Definition: DAkOmegaSSTFIML.H:74
Foam::DAkOmegaSSTFIML::epsilonByk
tmp< volScalarField::Internal > epsilonByk(const volScalarField &F1, const volTensorField &gradU) const
Definition: DAkOmegaSSTFIML.C:420
Foam::DAkOmegaSSTFIML::correctModelStates
virtual void correctModelStates(wordList &modelStates) const
update the turbulence state for DAStateInfo::regStates_
Definition: DAkOmegaSSTFIML.C:458
Foam::DAkOmegaSSTFIML::updateIntermediateVariables
virtual void updateIntermediateVariables()
update any intermediate variables that are dependent on state variables and are used in calcResiduals
Definition: DAkOmegaSSTFIML.C:605
Foam::DAkOmegaSSTFIML::DomegaEff
tmp< volScalarField > DomegaEff(const volScalarField &F1) const
Definition: DAkOmegaSSTFIML.H:142
U
U
Definition: pEqnPimpleDyM.H:60
Foam::DAkOmegaSSTFIML::tauRatio_
volScalarField tauRatio_
Definition: DAkOmegaSSTFIML.H:199
Foam::DAkOmegaSSTFIML::correctBoundaryConditions
virtual void correctBoundaryConditions()
update turbulence variable boundary values
Definition: DAkOmegaSSTFIML.C:516
Foam::DAkOmegaSSTFIML::b1_
dimensionedScalar b1_
Definition: DAkOmegaSSTFIML.H:71
Foam::DAUtility::pyCalcBeta
static void * pyCalcBeta
define a function pointer template for Python call back
Definition: DAUtility.H:129
forAll
forAll(pseudoP.boundaryField(), patchI)
Definition: solvePseudoPEqn.H:10
Foam::DATurbulenceModel::phi_
surfaceScalarField & phi_
face flux
Definition: DATurbulenceModel.H:89
Foam::DAUtility::pyCalcBetaInterface
static pyComputeInterface pyCalcBetaInterface
Definition: DAUtility.H:130
Foam::DAkOmegaSSTFIML::c1_
dimensionedScalar c1_
Definition: DAkOmegaSSTFIML.H:72
Foam::DAkOmegaSSTFIML::omegaSource
tmp< fvScalarMatrix > omegaSource() const
Definition: DAkOmegaSSTFIML.C:436
Foam::DAkOmegaSSTFIML::convectionTKE_
volScalarField convectionTKE_
Definition: DAkOmegaSSTFIML.H:198
Foam::DAkOmegaSSTFIML::omega_
volScalarField & omega_
Definition: DAkOmegaSSTFIML.H:177
Foam::DAkOmegaSSTFIML::F3
tmp< volScalarField > F3() const
Definition: DAkOmegaSSTFIML.C:381
Foam::DAkOmegaSSTFIML::Pk
tmp< volScalarField::Internal > Pk(const volScalarField::Internal &G) const
Definition: DAkOmegaSSTFIML.C:414
Foam::DAOption
Definition: DAOption.H:29
Foam::DAkOmegaSSTFIML::curvature_
volScalarField curvature_
Definition: DAkOmegaSSTFIML.H:201
daOption
DAOption daOption(mesh, pyOptions_)
Foam::DAkOmegaSSTFIML::correctStateResidualModelCon
virtual void correctStateResidualModelCon(List< List< word >> &stateCon) const
update the original variable connectivity for the adjoint state residuals in stateCon
Definition: DAkOmegaSSTFIML.C:616
Foam::DAkOmegaSSTFIML::outputs_
scalar * outputs_
Definition: DAkOmegaSSTFIML.H:79
normalizeResiduals
#define normalizeResiduals(resName)
Definition: DAMacroFunctions.H:28
Foam::DAkOmegaSSTFIML::addModelResidualCon
virtual void addModelResidualCon(HashTable< List< List< word >>> &allCon) const
add the model residual connectivity to stateCon
Definition: DAkOmegaSSTFIML.C:715
Foam::DAkOmegaSSTFIML::UGradMisalignment_
volScalarField UGradMisalignment_
Definition: DAkOmegaSSTFIML.H:202
Foam::DAkOmegaSSTFIML::GbyNu
tmp< volScalarField::Internal > GbyNu(const volScalarField::Internal &GbyNu0, const volScalarField::Internal &F2, const volScalarField::Internal &S2) const
Definition: DAkOmegaSSTFIML.C:403
Foam::DATurbulenceModel::phaseRhoPhi_
surfaceScalarField & phaseRhoPhi_
phase*phi*density field
Definition: DATurbulenceModel.H:95
Foam::DATurbulenceModel::allOptions_
const dictionary & allOptions_
all DAFoam options
Definition: DATurbulenceModel.H:80
Foam::DAkOmegaSSTFIML::a1_
dimensionedScalar a1_
Definition: DAkOmegaSSTFIML.H:70
Foam::DAkOmegaSSTFIML::betaStar_
dimensionedScalar betaStar_
Definition: DAkOmegaSSTFIML.H:68
Foam::DAkOmegaSSTFIML::pGradAlongStream_
volScalarField pGradAlongStream_
Definition: DAkOmegaSSTFIML.H:194
Foam::DAkOmegaSSTFIML::correct
virtual void correct()
solve the residual equations and update the state
Definition: DAkOmegaSSTFIML.C:807
Foam::DAkOmegaSSTFIML::Qsas
tmp< fvScalarMatrix > Qsas(const volScalarField::Internal &S2, const volScalarField::Internal &gamma, const volScalarField::Internal &beta) const
Definition: DAkOmegaSSTFIML.C:445
Foam::DAkOmegaSSTFIML::F2
tmp< volScalarField > F2() const
Definition: DAkOmegaSSTFIML.C:369
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAkOmegaSSTFIML::betaFieldInversionML_
volScalarField betaFieldInversionML_
Definition: DAkOmegaSSTFIML.H:191
Foam::DAkOmegaSSTFIML::calcResiduals
virtual void calcResiduals(const dictionary &options)
compute the turbulence residuals
Definition: DAkOmegaSSTFIML.C:972
Foam::DAkOmegaSSTFIML::QCriterion_
volScalarField QCriterion_
Definition: DAkOmegaSSTFIML.H:192
Foam::DAkOmegaSSTFIML::printInterval_
label printInterval_
time step interval to print residual
Definition: DAkOmegaSSTFIML.H:226
Foam::DAUtility::boundVar
static void boundVar(const dictionary &allOptions, volScalarField &var, const label printToScreen)
bound a volScalar variable based on parametes defined in DAOption::allOptions_
Definition: DAUtility.C:419
Foam::DAkOmegaSSTFIML::DkEff
tmp< volScalarField > DkEff(const volScalarField &F1) const
Definition: DAkOmegaSSTFIML.H:135
Foam::DAkOmegaSSTFIML::ReT_
volScalarField ReT_
Definition: DAkOmegaSSTFIML.H:197
Foam::DAkOmegaSSTFIML::F23
tmp< volScalarField > F23() const
Definition: DAkOmegaSSTFIML.C:391
Foam::DAkOmegaSSTFIML::calcBetaField
void calcBetaField()
calculate the beta field using the trained model
Definition: DAkOmegaSSTFIML.C:827
Foam::DATurbulenceModel::U_
volVectorField & U_
velocity
Definition: DATurbulenceModel.H:86
Foam::DATurbulenceModel::nut_
volScalarField & nut_
turbulence viscosity
Definition: DATurbulenceModel.H:83
Foam::DAkOmegaSSTFIML::saveOmegaNearWall
void saveOmegaNearWall()
save near wall omega values to omegaNearWall_
Definition: DAkOmegaSSTFIML.C:555
Foam::DAkOmegaSSTFIML::correctNut
virtual void correctNut()
update nut based on other turbulence variables and update the BCs
Definition: DAkOmegaSSTFIML.C:494
Foam::DAkOmegaSSTFIML::alphaOmega2_
dimensionedScalar alphaOmega2_
Definition: DAkOmegaSSTFIML.H:60
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DAkOmegaSSTFIML::inputs_
scalar * inputs_
inputs and outputs for the beta calculation
Definition: DAkOmegaSSTFIML.H:78
Foam::DAkOmegaSSTFIML::omegaNearWall_
scalarList omegaNearWall_
Definition: DAkOmegaSSTFIML.H:220
Foam::DAkOmegaSSTFIML::turbulenceIntensity_
volScalarField turbulenceIntensity_
Definition: DAkOmegaSSTFIML.H:195
Foam::DATurbulenceModel::nu
tmp< volScalarField > nu() const
get the nu field
Definition: DATurbulenceModel.C:262
Foam::DAkOmegaSSTFIML::setOmegaNearWall
void setOmegaNearWall()
set omegaNearWall_ to near wall omega values
Definition: DAkOmegaSSTFIML.C:579
Foam::DAkOmegaSSTFIML::kRes_
volScalarField kRes_
Definition: DAkOmegaSSTFIML.H:183
Foam::DATurbulenceModel::correctAlphat
void correctAlphat()
update alphat
Definition: DATurbulenceModel.C:206
Foam::DAkOmegaSSTFIML::k_
volScalarField & k_
Definition: DAkOmegaSSTFIML.H:182
solve
pseudoPEqn solve(solverDictP_)
Foam::DATurbulenceModel
Definition: DATurbulenceModel.H:61
Foam::DAkOmegaSSTFIML::gamma
tmp< volScalarField::Internal > gamma(const volScalarField::Internal &F1) const
Definition: DAkOmegaSSTFIML.H:128
Foam::DAkOmegaSSTFIML::pressureStress_
volScalarField pressureStress_
Definition: DAkOmegaSSTFIML.H:200
Foam::DAkOmegaSSTFIML::F1
tmp< volScalarField > F1(const volScalarField &CDkOmega) const
Definition: DAkOmegaSSTFIML.C:350
Foam::DAkOmegaSSTFIML::y_
const volScalarField & y_
Weight of betaFieldInversion for the turbulent transport equations if both are modified.
Definition: DAkOmegaSSTFIML.H:214
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFvSource, 0)
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFvSource, DAFvSourceActuatorDisk, dictionary)
Foam::DATurbulenceModel::mesh_
const fvMesh & mesh_
fvMesh
Definition: DATurbulenceModel.H:74
Foam::DAkOmegaSSTFIML::p_
volScalarField & p_
Definition: DAkOmegaSSTFIML.H:193
Foam::DAkOmegaSSTFIML::solveTurbState_
label solveTurbState_
whether to solve for turb states
Definition: DAkOmegaSSTFIML.H:223
Foam::DATurbulenceModel::isPrintTime
label isPrintTime(const Time &runTime, const label printInterval) const
Definition: DATurbulenceModel.C:464
DAkOmegaSSTFIML.H
Foam::DAkOmegaSSTFIML::omegaRes_
volScalarField omegaRes_
Definition: DAkOmegaSSTFIML.H:178
Foam::DAkOmegaSSTFIML::beta
tmp< volScalarField::Internal > beta(const volScalarField::Internal &F1) const
Definition: DAkOmegaSSTFIML.H:122
Foam::DAkOmegaSSTFIML::kSource
tmp< fvScalarMatrix > kSource() const
Definition: DAkOmegaSSTFIML.C:427
Foam::DATurbulenceModel::phase_
volScalarField phase_
phase field
Definition: DATurbulenceModel.H:92
Foam::DAkOmegaSSTFIML::correctOmegaBoundaryConditions
void correctOmegaBoundaryConditions()
specially treatment to correct epsilon BC
Definition: DAkOmegaSSTFIML.C:528
rho
volScalarField & rho
Definition: createRefsRhoSimpleC.H:8