DAkOmega.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 "DAkOmega.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 defineTypeNameAndDebug(DAkOmega, 0);
38 addToRunTimeSelectionTable(DATurbulenceModel, DAkOmega, dictionary);
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
42  const word modelType,
43  const fvMesh& mesh,
44  const DAOption& daOption)
45  : DATurbulenceModel(modelType, mesh, daOption),
46  // kOmega parameters
47  Cmu_(dimensioned<scalar>::lookupOrAddToDict(
48  "betaStar",
49  this->coeffDict_,
50  0.09)),
51  beta_(dimensioned<scalar>::lookupOrAddToDict(
52  "beta",
53  this->coeffDict_,
54  0.072)),
55  gamma_(dimensioned<scalar>::lookupOrAddToDict(
56  "gamma",
57  this->coeffDict_,
58  0.52)),
59  alphaK_(dimensioned<scalar>::lookupOrAddToDict(
60  "alphaK",
61  this->coeffDict_,
62  0.5)),
63  alphaOmega_(dimensioned<scalar>::lookupOrAddToDict(
64  "alphaOmega",
65  this->coeffDict_,
66  0.5)),
67  // Augmented variables
68  omega_(const_cast<volScalarField&>(
69  mesh_.thisDb().lookupObject<volScalarField>("omega"))),
70  omegaRes_(
71  IOobject(
72  "omegaRes",
73  mesh.time().timeName(),
74  mesh,
75  IOobject::NO_READ,
76  IOobject::NO_WRITE),
77  mesh,
78 #ifdef CompressibleFlow
79  dimensionedScalar("omegaRes", dimensionSet(1, -3, -2, 0, 0, 0, 0), 0.0),
80 #endif
81 #ifdef IncompressibleFlow
82  dimensionedScalar("omegaRes", dimensionSet(0, 0, -2, 0, 0, 0, 0), 0.0),
83 #endif
84  zeroGradientFvPatchField<scalar>::typeName),
85  k_(
86  const_cast<volScalarField&>(
87  mesh_.thisDb().lookupObject<volScalarField>("k"))),
88  kRes_(
89  IOobject(
90  "kRes",
91  mesh.time().timeName(),
92  mesh,
93  IOobject::NO_READ,
94  IOobject::NO_WRITE),
95  mesh,
96 #ifdef CompressibleFlow
97  dimensionedScalar("kRes", dimensionSet(1, -1, -3, 0, 0, 0, 0), 0.0),
98 #endif
99 #ifdef IncompressibleFlow
100  dimensionedScalar("kRes", dimensionSet(0, 2, -3, 0, 0, 0, 0), 0.0),
101 #endif
102  zeroGradientFvPatchField<scalar>::typeName)
103 {
104 
105  // initialize printInterval_ we need to check whether it is a steady state
106  // or unsteady primal solver
107  IOdictionary fvSchemes(
108  IOobject(
109  "fvSchemes",
110  mesh.time().system(),
111  mesh,
112  IOobject::MUST_READ,
113  IOobject::NO_WRITE,
114  false));
115  word ddtScheme = word(fvSchemes.subDict("ddtSchemes").lookup("default"));
116  if (ddtScheme == "steadyState")
117  {
119  daOption.getAllOptions().lookupOrDefault<label>("printInterval", 100);
120  }
121  else
122  {
124  daOption.getAllOptions().lookupOrDefault<label>("printIntervalUnsteady", 500);
125  }
126 
127  // calculate the size of omegaWallFunction faces
128  label nWallFaces = 0;
129  forAll(omega_.boundaryField(), patchI)
130  {
131  if (omega_.boundaryField()[patchI].type() == "omegaWallFunction"
132  && omega_.boundaryField()[patchI].size() > 0)
133  {
134  forAll(omega_.boundaryField()[patchI], faceI)
135  {
136  nWallFaces++;
137  }
138  }
139  }
140 
141  // initialize omegaNearWall
142  omegaNearWall_.setSize(nWallFaces);
143 }
144 
145 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
146 
147 // Augmented functions
148 void DAkOmega::correctModelStates(wordList& modelStates) const
149 {
150  /*
151  Description:
152  Update the name in modelStates based on the selected physical model at runtime
153 
154  Example:
155  In DAStateInfo, if the modelStates reads:
156 
157  modelStates = {"nut"}
158 
159  then for the SA model, calling correctModelStates(modelStates) will give:
160 
161  modelStates={"nuTilda"}
162 
163  while calling correctModelStates(modelStates) for the SST model will give
164 
165  modelStates={"k","omega"}
166 
167  We don't udpate the names for the radiation model because users are
168  supposed to set modelStates={"G"}
169  */
170 
171  // For SST model, we need to replace nut with k, omega
172 
173  forAll(modelStates, idxI)
174  {
175  word stateName = modelStates[idxI];
176  if (stateName == "nut")
177  {
178  modelStates[idxI] = "omega";
179  modelStates.append("k");
180  }
181  }
182 }
183 
185 {
186  /*
187  Description:
188  Update nut based on other turbulence variables and update the BCs
189  Also update alphat if is present
190  */
191 
192  nut_ = k_ / omega_;
193 
194  nut_.correctBoundaryConditions(); // nutkWallFunction: update wall face nut based on k
195 
196  // this is basically BasicTurbulenceModel::correctNut();
197  this->correctAlphat();
198 
199  return;
200 }
201 
203 {
204  /*
205  Description:
206  Update turbulence variable boundary values
207  */
208 
209  // correct the BCs for the perturbed fields
210  // kqWallFunction is a zero-gradient BC
211  k_.correctBoundaryConditions();
212 }
213 
215 {
216  /*
217  Description:
218  this is a special treatment for omega BC because we cant directly call omega.
219  correctBoundaryConditions() because it will modify the internal omega and G that
220  are right next to walls. This will mess up adjoint Jacobians
221  To solve this issue,
222  1. we store the near wall omega before calling omega.correctBoundaryConditions()
223  2. call omega.correctBoundaryConditions()
224  3. Assign the stored near wall omega back
225  4. Apply a zeroGradient BC for omega at the wall patches
226  *********** NOTE *************
227  this treatment will obviously downgrade the accuracy of adjoint derivative since it is
228  not 100% consistent with what is used for the flow solver; however, our observation
229  shows that the impact is not very large.
230  */
231 
232  // save the perturbed omega at the wall
233  this->saveOmegaNearWall();
234  // correct omega boundary conditions, this includes updating wall face and near wall omega values,
235  // updating the inter-proc BCs
236  omega_.correctBoundaryConditions();
237  // reset the corrected omega near wall cell to its perturbed value
238  this->setOmegaNearWall();
239 }
240 
242 {
243  /*
244  Description:
245  Save the current near wall omega values to omegaNearWall_
246  */
247  label counterI = 0;
248  forAll(omega_.boundaryField(), patchI)
249  {
250  if (omega_.boundaryField()[patchI].type() == "omegaWallFunction"
251  and omega_.boundaryField()[patchI].size() > 0)
252  {
253  const UList<label>& faceCells = mesh_.boundaryMesh()[patchI].faceCells();
254  forAll(faceCells, faceI)
255  {
256  //Info<<"Near Wall cellI: "<<faceCells[faceI]<<endl;
257  omegaNearWall_[counterI] = omega_[faceCells[faceI]];
258  counterI++;
259  }
260  }
261  }
262  return;
263 }
264 
266 {
267  /*
268  Description:
269  Set the current near wall omega values from omegaNearWall_
270  Here we also apply a zeroGradient BC to the wall faces
271  */
272  label counterI = 0;
273  forAll(omega_.boundaryField(), patchI)
274  {
275  if (omega_.boundaryField()[patchI].type() == "omegaWallFunction"
276  && omega_.boundaryField()[patchI].size() > 0)
277  {
278  const UList<label>& faceCells = mesh_.boundaryMesh()[patchI].faceCells();
279  forAll(faceCells, faceI)
280  {
281  omega_[faceCells[faceI]] = omegaNearWall_[counterI];
282  // zeroGradient BC
283  omega_.boundaryFieldRef()[patchI][faceI] = omega_[faceCells[faceI]];
284  counterI++;
285  }
286  }
287  }
288  return;
289 }
290 
292 {
293  /*
294  Description:
295  Update nut based on nuTilda. Note: we need to update nut and its BC since we
296  may have perturbed other turbulence vars that affect the nut values
297  */
298 
299  this->correctNut();
300 }
301 
302 void DAkOmega::correctStateResidualModelCon(List<List<word>>& stateCon) const
303 {
304  /*
305  Description:
306  Update the original variable connectivity for the adjoint state
307  residuals in stateCon. Basically, we modify/add state variables based on the
308  original model variables defined in stateCon.
309 
310  Input:
311 
312  stateResCon: the connectivity levels for a state residual, defined in Foam::DAJacCon
313 
314  Example:
315  If stateCon reads:
316  stateCon=
317  {
318  {"U", "p", "nut"},
319  {"p"}
320  }
321 
322  For the SA turbulence model, calling this function for will get a new stateCon
323  stateCon=
324  {
325  {"U", "p", "nuTilda"},
326  {"p"}
327  }
328 
329  For the SST turbulence model, calling this function will give
330  stateCon=
331  {
332  {"U", "p", "k", "omega"},
333  {"p", "U"}
334  }
335  ***NOTE***: we add a extra level of U connectivity because nut is
336  related to grad(U), k, and omega in SST!
337  */
338 
339  label stateConSize = stateCon.size();
340  forAll(stateCon, idxI)
341  {
342  label addUCon = 0;
343  forAll(stateCon[idxI], idxJ)
344  {
345  word conStateName = stateCon[idxI][idxJ];
346  if (conStateName == "nut")
347  {
348  stateCon[idxI][idxJ] = "omega"; // replace nut with omega
349  stateCon[idxI].append("k"); // also add k for that level
350  addUCon = 1;
351  }
352  }
353  // add U for the current level and level+1 if it is not there yet
354  label isU;
355  if (addUCon == 1)
356  {
357  // first add U for the current level
358  isU = 0;
359  forAll(stateCon[idxI], idxJ)
360  {
361  word conStateName = stateCon[idxI][idxJ];
362  if (conStateName == "U")
363  {
364  isU = 1;
365  }
366  }
367  if (!isU)
368  {
369  stateCon[idxI].append("U");
370  }
371 
372  // now add U for level+1 if idxI is not the largest level
373  // if idxI is already the largest level, we have a problem
374  if (idxI != stateConSize - 1)
375  {
376  isU = 0;
377  forAll(stateCon[idxI + 1], idxJ)
378  {
379  word conStateName = stateCon[idxI + 1][idxJ];
380  if (conStateName == "U")
381  {
382  isU = 1;
383  }
384  }
385  if (!isU)
386  {
387  stateCon[idxI + 1].append("U");
388  }
389  }
390  else
391  {
392  FatalErrorIn(
393  "In DAStateInfo, nut shows in the largest connectivity level! "
394  "This is not supported!")
395  << exit(FatalError);
396  }
397  }
398  }
399 }
400 
401 void DAkOmega::addModelResidualCon(HashTable<List<List<word>>>& allCon) const
402 {
403  /*
404  Description:
405  Add the connectivity levels for all physical model residuals to allCon
406 
407  Input:
408  allCon: the connectivity levels for all state residual, defined in DAJacCon
409 
410  Example:
411  If stateCon reads:
412  allCon=
413  {
414  "URes":
415  {
416  {"U", "p", "nut"},
417  {"p"}
418  }
419  }
420 
421  For the SA turbulence model, calling this function for will get a new stateCon,
422  something like this:
423  allCon=
424  {
425  "URes":
426  {
427  {"U", "p", "nuTilda"},
428  {"p"}
429  },
430  "nuTildaRes":
431  {
432  {"U", "phi", "nuTilda"},
433  {"U"}
434  }
435  }
436 
437  */
438 
439  word pName;
440 
441  if (mesh_.thisDb().foundObject<volScalarField>("p"))
442  {
443  pName = "p";
444  }
445  else if (mesh_.thisDb().foundObject<volScalarField>("p_rgh"))
446  {
447  pName = "p_rgh";
448  }
449  else
450  {
451  FatalErrorIn(
452  "Neither p nor p_rgh was found in mesh.thisDb()!"
453  "addModelResidualCon failed to setup turbulence residuals!")
454  << exit(FatalError);
455  }
456 
457  // NOTE: for compressible flow, it depends on rho so we need to add T and p
458 #ifdef IncompressibleFlow
459  allCon.set(
460  "omegaRes",
461  {
462  {"U", "omega", "k", "phi"}, // lv0
463  {"U", "omega", "k"}, // lv1
464  {"U", "omega", "k"} // lv2
465  });
466  allCon.set(
467  "kRes",
468  {
469  {"U", "omega", "k", "phi"}, // lv0
470  {"U", "omega", "k"}, // lv1
471  {"U", "omega", "k"} // lv2
472  });
473 #endif
474 
475 #ifdef CompressibleFlow
476  allCon.set(
477  "omegaRes",
478  {
479  {"U", "T", pName, "omega", "k", "phi"}, // lv0
480  {"U", "T", pName, "omega", "k"}, // lv1
481  {"U", "T", pName, "omega", "k"} // lv2
482  });
483  allCon.set(
484  "kRes",
485  {
486  {"U", "T", pName, "omega", "k", "phi"}, // lv0
487  {"U", "T", pName, "omega", "k"}, // lv1
488  {"U", "T", pName, "omega", "k"} // lv2
489  });
490 #endif
491 }
492 
494 {
495  /*
496  Descroption:
497  Solve the residual equations and update the state. This function will be called
498  by the DASolver. It is needed because we want to control the output frequency
499  of the residual convergence every 100 steps. If using the correct from turbulence
500  it will output residual every step which will be too much of information.
501  */
502 
503  // We set the flag solveTurbState_ to 1 such that in the calcResiduals function
504  // we will solve and update nuTilda
505  solveTurbState_ = 1;
506  dictionary dummyOptions;
507  this->calcResiduals(dummyOptions);
508  // after it, we reset solveTurbState_ = 0 such that calcResiduals will not
509  // update nuTilda when calling from the adjoint class, i.e., solveAdjoint from DASolver.
510  solveTurbState_ = 0;
511 }
512 
513 void DAkOmega::calcResiduals(const dictionary& options)
514 {
515  /*
516  Descroption:
517  If solveTurbState_ == 1, this function solve and update k and omega, and
518  is the same as calling turbulence.correct(). If solveTurbState_ == 0,
519  this function compute residuals for turbulence variables, e.g., nuTildaRes_
520 
521  Input:
522  options.isPC: 1 means computing residuals for preconditioner matrix.
523  This essentially use the first order scheme for div(phi,nuTilda)
524 
525  p_, U_, phi_, etc: State variables in OpenFOAM
526 
527  Output:
528  kRes_/omegaRes_: If solveTurbState_ == 0, update the residual field variable
529 
530  k_/omega_: If solveTurbState_ == 1, update them
531  */
532 
533  // Copy and modify based on the "correct" function
534 
535  label printToScreen = this->isPrintTime(mesh_.time(), printInterval_);
536 
537  word divKScheme = "div(phi,k)";
538  word divOmegaScheme = "div(phi,omega)";
539 
540  label isPC = 0;
541 
542  if (!solveTurbState_)
543  {
544  isPC = options.getLabel("isPC");
545 
546  if (isPC)
547  {
548  divKScheme = "div(pc)";
549  divOmegaScheme = "div(pc)";
550  }
551  }
552 
553  // Note: for compressible flow, the "this->phi()" function divides phi by fvc:interpolate(rho),
554  // while for the incompresssible "this->phi()" returns phi only
555  // see src/TurbulenceModels/compressible/compressibleTurbulenceModel.C line 62 to 73
556  volScalarField divU(fvc::div(fvc::absolute(phi_ / fvc::interpolate(rho_), U_)));
557 
558  tmp<volTensorField> tgradU = fvc::grad(U_);
559  volScalarField G("kOmega:G", nut_ * (tgradU() && dev(twoSymm(tgradU()))));
560  tgradU.clear();
561 
562  if (solveTurbState_)
563  {
564  // Update omega and G at the wall
565  omega_.boundaryFieldRef().updateCoeffs();
566  }
567  else
568  {
569  // NOTE instead of calling omega_.boundaryFieldRef().updateCoeffs();
570  // here we call our self-defined boundary conditions
572  }
573 
574  // Turbulent frequency equation
575  tmp<fvScalarMatrix> omegaEqn(
576  fvm::ddt(phase_, rho_, omega_)
577  + fvm::div(phaseRhoPhi_, omega_, divOmegaScheme)
578  - fvm::laplacian(phase_ * rho_ * DomegaEff(), omega_)
579  == gamma_ * phase_ * rho_ * G * omega_ / k_
580  - fvm::SuSp(scalar(2.0 / 3.0) * gamma_ * phase_ * rho_ * divU, omega_)
581  - fvm::Sp(beta_ * phase_ * rho_ * omega_, omega_));
582 
583  omegaEqn.ref().relax();
584  omegaEqn.ref().boundaryManipulate(omega_.boundaryFieldRef());
585 
586  if (solveTurbState_)
587  {
588 
589  // get the solver performance info such as initial
590  // and final residuals
591  SolverPerformance<scalar> solverOmega = solve(omegaEqn);
592  if (printToScreen)
593  {
594  Info << "omega Initial residual: " << solverOmega.initialResidual() << endl
595  << " Final residual: " << solverOmega.finalResidual() << endl;
596  }
597 
598  DAUtility::boundVar(allOptions_, omega_, printToScreen);
599  }
600  else
601  {
602  // reset the corrected omega near wall cell to its perturbed value
603  this->setOmegaNearWall();
604 
605  // calculate residuals
606  omegaRes_ = omegaEqn() & omega_;
607  // need to normalize residuals
608  normalizeResiduals(omegaRes);
609  }
610 
611  // Turbulent kinetic energy equation
612  tmp<fvScalarMatrix> kEqn(
613  fvm::ddt(phase_, rho_, k_)
614  + fvm::div(phaseRhoPhi_, k_, divKScheme)
615  - fvm::laplacian(phase_ * rho_ * DkEff(), k_)
616  == phase_ * rho_ * G
617  - fvm::SuSp((2.0 / 3.0) * phase_ * rho_ * divU, k_)
618  - fvm::Sp(Cmu_ * phase_ * rho_ * omega_, k_));
619 
620  kEqn.ref().relax();
621 
622  if (solveTurbState_)
623  {
624 
625  // get the solver performance info such as initial
626  // and final residuals
627  SolverPerformance<scalar> solverK = solve(kEqn);
628  if (printToScreen)
629  {
630  Info << "k Initial residual: " << solverK.initialResidual() << endl
631  << " Final residual: " << solverK.finalResidual() << endl;
632  }
633 
634  DAUtility::boundVar(allOptions_, k_, printToScreen);
635 
636  this->correctNut();
637  }
638  else
639  {
640  // calculate residuals
641  kRes_ = kEqn() & k_;
642  // need to normalize residuals
643  normalizeResiduals(kRes);
644  }
645 
646  return;
647 }
648 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
649 
650 } // End namespace Foam
651 
652 // ************************************************************************* //
forAll
forAll(pseudoP.boundaryField(), patchI)
Definition: solvePseudoPEqn.H:10
Foam::DATurbulenceModel::phi_
surfaceScalarField & phi_
face flux
Definition: DATurbulenceModel.H:89
Foam::DAkOmega::saveOmegaNearWall
void saveOmegaNearWall()
save near wall omega values to omegaNearWall_
Definition: DAkOmega.C:241
Foam::DAOption
Definition: DAOption.H:29
daOption
DAOption daOption(mesh, pyOptions_)
Foam::DAkOmega::Cmu_
dimensionedScalar Cmu_
Definition: DAkOmega.H:55
normalizeResiduals
#define normalizeResiduals(resName)
Definition: DAMacroFunctions.H:28
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::DAkOmega::omegaRes_
volScalarField omegaRes_
Definition: DAkOmega.H:94
Foam::DAkOmega::correct
virtual void correct()
solve the residual equations and update the state
Definition: DAkOmega.C:493
Foam::DAkOmega::DomegaEff
tmp< volScalarField > DomegaEff() const
Definition: DAkOmega.H:78
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAkOmega::kRes_
volScalarField kRes_
Definition: DAkOmega.H:96
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::DAkOmega::omegaNearWall_
scalarList omegaNearWall_
Definition: DAkOmega.H:103
Foam::DATurbulenceModel::U_
volVectorField & U_
velocity
Definition: DATurbulenceModel.H:86
Foam::DATurbulenceModel::nut_
volScalarField & nut_
turbulence viscosity
Definition: DATurbulenceModel.H:83
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DAkOmega::correctStateResidualModelCon
virtual void correctStateResidualModelCon(List< List< word >> &stateCon) const
update the original variable connectivity for the adjoint state residuals in stateCon
Definition: DAkOmega.C:302
Foam::DAkOmega::k_
volScalarField & k_
Definition: DAkOmega.H:95
Foam::DAkOmega::correctModelStates
virtual void correctModelStates(wordList &modelStates) const
update the turbulence state for DAStateInfo::regStates_
Definition: DAkOmega.C:148
Foam::DAkOmega::solveTurbState_
label solveTurbState_
whether to solve for turb states
Definition: DAkOmega.H:106
Foam::DAkOmega::omega_
volScalarField & omega_
Definition: DAkOmega.H:93
Foam::DAkOmega::gamma_
dimensionedScalar gamma_
Definition: DAkOmega.H:57
Foam::DAkOmega::correctNut
virtual void correctNut()
update nut based on other turbulence variables and update the BCs
Definition: DAkOmega.C:184
Foam::DATurbulenceModel::correctAlphat
void correctAlphat()
update alphat
Definition: DATurbulenceModel.C:206
Foam::DAkOmega::addModelResidualCon
virtual void addModelResidualCon(HashTable< List< List< word >>> &allCon) const
add the model residual connectivity to stateCon
Definition: DAkOmega.C:401
solve
pseudoPEqn solve(solverDictP_)
Foam::DATurbulenceModel
Definition: DATurbulenceModel.H:61
Foam::DAkOmega::beta_
dimensionedScalar beta_
Definition: DAkOmega.H:56
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFvSource, 0)
Foam::DAkOmega::calcResiduals
virtual void calcResiduals(const dictionary &options)
compute the turbulence residuals
Definition: DAkOmega.C:513
Foam::DAkOmega::correctOmegaBoundaryConditions
void correctOmegaBoundaryConditions()
specially treatment to correct epsilon BC
Definition: DAkOmega.C:214
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFvSource, DAFvSourceActuatorDisk, dictionary)
Foam::DATurbulenceModel::mesh_
const fvMesh & mesh_
fvMesh
Definition: DATurbulenceModel.H:74
Foam::DATurbulenceModel::isPrintTime
label isPrintTime(const Time &runTime, const label printInterval) const
Definition: DATurbulenceModel.C:464
DAkOmega.H
Foam::DAkOmega::printInterval_
label printInterval_
time step interval to print residual
Definition: DAkOmega.H:109
Foam::DAkOmega::correctBoundaryConditions
virtual void correctBoundaryConditions()
update turbulence variable boundary values
Definition: DAkOmega.C:202
Foam::DAkOmega::updateIntermediateVariables
virtual void updateIntermediateVariables()
update any intermediate variables that are dependent on state variables and are used in calcResiduals
Definition: DAkOmega.C:291
Foam::DAkOmega::setOmegaNearWall
void setOmegaNearWall()
set omegaNearWall_ to near wall omega values
Definition: DAkOmega.C:265
Foam::DATurbulenceModel::phase_
volScalarField phase_
phase field
Definition: DATurbulenceModel.H:92
Foam::DAkOmega::DkEff
tmp< volScalarField > DkEff() const
Definition: DAkOmega.H:65
Foam::DAkOmega::DAkOmega
DAkOmega(const word modelType, const fvMesh &mesh, const DAOption &daOption)
Definition: DAkOmega.C:41