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