DASpalartAllmaras.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v4
5 
6  This file is modified from OpenFOAM's source code
7  src/TurbulenceModels/turbulenceModels/RAS/SpalartAllmaras/SpalartAllmaras.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 "DASpalartAllmaras.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 defineTypeNameAndDebug(DASpalartAllmaras, 0);
38 addToRunTimeSelectionTable(DATurbulenceModel, DASpalartAllmaras, dictionary);
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
42  const word modelType,
43  const fvMesh& mesh,
44  const DAOption& daOption)
45  : DATurbulenceModel(modelType, mesh, daOption),
46  // SA parameters
47  sigmaNut_(dimensioned<scalar>::lookupOrAddToDict(
48  "sigmaNut",
49  this->coeffDict_,
50  0.66666)),
51  kappa_(dimensioned<scalar>::lookupOrAddToDict(
52  "kappa",
53  this->coeffDict_,
54  0.41)),
55 
56  Cb1_(dimensioned<scalar>::lookupOrAddToDict(
57  "Cb1",
58  this->coeffDict_,
59  0.1355)),
60  Cb2_(dimensioned<scalar>::lookupOrAddToDict(
61  "Cb2",
62  this->coeffDict_,
63  0.622)),
64  Cw1_(Cb1_ / sqr(kappa_) + (1.0 + Cb2_) / sigmaNut_),
65  Cw2_(dimensioned<scalar>::lookupOrAddToDict(
66  "Cw2",
67  this->coeffDict_,
68  0.3)),
69  Cw3_(dimensioned<scalar>::lookupOrAddToDict(
70  "Cw3",
71  this->coeffDict_,
72  2.0)),
73  Cv1_(dimensioned<scalar>::lookupOrAddToDict(
74  "Cv1",
75  this->coeffDict_,
76  7.1)),
77  Cs_(dimensioned<scalar>::lookupOrAddToDict(
78  "Cs",
79  this->coeffDict_,
80  0.3)),
81  // Augmented variables
82  nuTilda_(const_cast<volScalarField&>(
83  mesh.thisDb().lookupObject<volScalarField>("nuTilda"))),
84  nuTildaRes_(
85  IOobject(
86  "nuTildaRes",
87  mesh.time().timeName(),
88  mesh,
89  IOobject::NO_READ,
90  IOobject::NO_WRITE),
91  mesh,
92  dimensionedScalar("nuTildaRes", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0),
93  zeroGradientFvPatchField<scalar>::typeName),
94  y_(mesh.thisDb().lookupObject<volScalarField>("yWall")),
95  betaFINuTilda_(
96  IOobject(
97  "betaFINuTilda",
98  mesh.time().timeName(),
99  mesh,
100  IOobject::READ_IF_PRESENT,
101  IOobject::AUTO_WRITE),
102  mesh,
103  dimensionedScalar("betaFINuTilda", dimensionSet(0, 0, 0, 0, 0, 0, 0), 1.0),
104  "zeroGradient"),
105  psiNuTildaPC_(nullptr),
106  dPsiNuTilda_("dPsiNuTilda", nuTilda_)
107 {
108  // we need to reset the nuTildaRes's dimension based on the model type
109  if (turbModelType_ == "incompressible")
110  {
111  nuTildaRes_.dimensions().reset(dimensionSet(0, 2, -2, 0, 0, 0, 0));
112  }
113 
114  if (turbModelType_ == "compressible")
115  {
116  nuTildaRes_.dimensions().reset(dimensionSet(1, -1, -2, 0, 0, 0, 0));
117  }
118 }
119 
120 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
121 
122 // SA member functions. these functions are copied from
123 // src/TurbulenceModels/turbulenceModels/RAS/SpalartAllmaras/SpalartAllmaras.C
124 tmp<volScalarField> DASpalartAllmaras::chi() const
125 {
126  return nuTilda_ / this->nu();
127 }
128 
129 tmp<volScalarField> DASpalartAllmaras::fv1(
130  const volScalarField& chi) const
131 {
132  const volScalarField chi3(pow3(chi));
133  return chi3 / (chi3 + pow3(Cv1_));
134 }
135 
136 tmp<volScalarField> DASpalartAllmaras::fv2(
137  const volScalarField& chi,
138  const volScalarField& fv1) const
139 {
140  return 1.0 - chi / (1.0 + chi * fv1);
141 }
142 
143 tmp<volScalarField> DASpalartAllmaras::Stilda(
144  const volScalarField& chi,
145  const volScalarField& fv1) const
146 {
147  volScalarField Omega(::sqrt(2.0) * mag(skew(fvc::grad(U_))));
148 
149  return (
150  max(
151  Omega
152  + fv2(chi, fv1) * nuTilda_ / sqr(kappa_ * y_),
153  Cs_ * Omega));
154 }
155 
156 tmp<volScalarField> DASpalartAllmaras::fw(
157  const volScalarField& Stilda) const
158 {
159  volScalarField r(
160  min(
161  nuTilda_
162  / (max(
163  Stilda,
164  dimensionedScalar("SMALL", Stilda.dimensions(), SMALL))
165  * sqr(kappa_ * y_)),
166  scalar(10.0)));
167  r.boundaryFieldRef() == 0.0;
168 
169  const volScalarField g(r + Cw2_ * (pow6(r) - r));
170 
171  return g * pow((1.0 + pow6(Cw3_)) / (pow6(g) + pow6(Cw3_)), 1.0 / 6.0);
172 }
173 
174 tmp<volScalarField> DASpalartAllmaras::DnuTildaEff() const
175 {
176  return tmp<volScalarField>(
177  new volScalarField("DnuTildaEff", (nuTilda_ + this->nu()) / sigmaNut_));
178 }
179 
180 // Augmented functions
181 void DASpalartAllmaras::correctModelStates(wordList& modelStates) const
182 {
183  /*
184  Description:
185  Update the name in modelStates based on the selected physical model at runtime
186 
187  Example:
188  In DAStateInfo, if the modelStates reads:
189 
190  modelStates = {"nut"}
191 
192  then for the SA model, calling correctModelStates(modelStates) will give:
193 
194  modelStates={"nuTilda"}
195 
196  while calling correctModelStates(modelStates) for the SST model will give
197 
198  modelStates={"k","omega"}
199 
200  We don't udpate the names for the radiation model becasue users are
201  supposed to set modelStates={"G"}
202  */
203 
204  // replace nut with nuTilda
205  forAll(modelStates, idxI)
206  {
207  word stateName = modelStates[idxI];
208  if (stateName == "nut")
209  {
210  modelStates[idxI] = "nuTilda";
211  }
212  }
213 }
214 
216 {
217  /*
218  Description:
219  Update nut based on other turbulence variables and update the BCs
220  Also update alphat if is present
221  */
222 
223  const volScalarField chi(this->chi());
224  const volScalarField fv1(this->fv1(chi));
225  nut_ = nuTilda_ * fv1;
226 
227  nut_.correctBoundaryConditions();
228 
229  // this is basically BasicTurbulenceModel::correctNut();
230  this->correctAlphat();
231 
232  return;
233 }
234 
236 {
237  /*
238  Description:
239  Update turbulence variable boundary values
240  */
241 
242  // correct the BCs for the perturbed fields
243  nuTilda_.correctBoundaryConditions();
244 }
245 
247 {
248  /*
249  Description:
250  Update nut based on nuTilda. Note: we need to update nut and its BC since we
251  may have perturbed other turbulence vars that affect the nut values
252  */
253 
254  this->correctNut();
255 }
256 
257 void DASpalartAllmaras::correctStateResidualModelCon(List<List<word>>& stateCon) const
258 {
259  /*
260  Description:
261  Update the original variable connectivity for the adjoint state
262  residuals in stateCon. Basically, we modify/add state variables based on the
263  original model variables defined in stateCon.
264 
265  Input:
266 
267  stateResCon: the connectivity levels for a state residual, defined in Foam::DAJacCon
268 
269  Example:
270  If stateCon reads:
271  stateCon=
272  {
273  {"U", "p", "nut"},
274  {"p"}
275  }
276 
277  For the SA turbulence model, calling this function for will get a new stateCon
278  stateCon=
279  {
280  {"U", "p", "nuTilda"},
281  {"p"}
282  }
283 
284  For the SST turbulence model, calling this function will give
285  stateCon=
286  {
287  {"U", "p", "k", "omega"},
288  {"p", "U"}
289  }
290  ***NOTE***: we add a extra level of U connectivity because nut is
291  related to grad(U), k, and omega in SST!
292  */
293 
294  forAll(stateCon, idxI)
295  {
296  forAll(stateCon[idxI], idxJ)
297  {
298  word conStateName = stateCon[idxI][idxJ];
299  if (conStateName == "nut")
300  {
301  stateCon[idxI][idxJ] = "nuTilda";
302  }
303  }
304  }
305 }
306 
307 void DASpalartAllmaras::addModelResidualCon(HashTable<List<List<word>>>& allCon) const
308 {
309  /*
310  Description:
311  Add the connectivity levels for all physical model residuals to allCon
312 
313  Input:
314  allCon: the connectivity levels for all state residual, defined in DAJacCon
315 
316  Example:
317  If stateCon reads:
318  allCon=
319  {
320  "URes":
321  {
322  {"U", "p", "nut"},
323  {"p"}
324  }
325  }
326 
327  For the SA turbulence model, calling this function for will get a new stateCon,
328  something like this:
329  allCon=
330  {
331  "URes":
332  {
333  {"U", "p", "nuTilda"},
334  {"p"}
335  },
336  "nuTildaRes":
337  {
338  {"U", "phi", "nuTilda"},
339  {"U"}
340  }
341  }
342 
343  */
344 
345  word pName;
346 
347  if (mesh_.thisDb().foundObject<volScalarField>("p"))
348  {
349  pName = "p";
350  }
351  else if (mesh_.thisDb().foundObject<volScalarField>("p_rgh"))
352  {
353  pName = "p_rgh";
354  }
355  else
356  {
357  FatalErrorIn(
358  "Neither p nor p_rgh was found in mesh.thisDb()!"
359  "addModelResidualCon failed to setup turbulence residuals!")
360  << exit(FatalError);
361  }
362 
363  // NOTE: for compressible flow, it depends on rho so we need to add T and p
364  if (turbModelType_ == "incompressible")
365  {
366  allCon.set(
367  "nuTildaRes",
368  {
369  {"U", "nuTilda", "phi"}, // lv0
370  {"U", "nuTilda"}, // lv1
371  {"nuTilda"} // lv2
372  });
373  }
374  else if (turbModelType_ == "compressible")
375  {
376  allCon.set(
377  "nuTildaRes",
378  {
379  {"U", "T", pName, "nuTilda", "phi"}, // lv0
380  {"U", "T", pName, "nuTilda"}, // lv1
381  {"T", pName, "nuTilda"} // lv2
382  });
383  }
384 }
385 
386 void DASpalartAllmaras::correct(label printToScreen)
387 {
388  /*
389  Descroption:
390  Solve the residual equations and update the state. This function will be called
391  by the DASolver. It is needed because we want to control the output frequency
392  of the residual convergence every 100 steps. If using the correct from turbulence
393  it will output residual every step which will be too much of information.
394  */
395 
396  // We set the flag solveTurbState_ to 1 such that in the calcResiduals function
397  // we will solve and update nuTilda
398  solveTurbState_ = 1;
399  dictionary dummyOptions;
400  dummyOptions.set("printToScreen", printToScreen);
401  this->calcResiduals(dummyOptions);
402  // after it, we reset solveTurbState_ = 0 such that calcResiduals will not
403  // update nuTilda when calling from the adjoint class, i.e., solveAdjoint from DASolver.
404  solveTurbState_ = 0;
405 }
406 
407 void DASpalartAllmaras::calcResiduals(const dictionary& options)
408 {
409  /*
410  Descroption:
411  If solveTurbState_ == 1, this function solve and update nuTilda, and
412  is the same as calling turbulence.correct(). If solveTurbState_ == 0,
413  this function compute residuals for turbulence variables, e.g., nuTildaRes_
414 
415  Input:
416  options.isPC: 1 means computing residuals for preconditioner matrix.
417  This essentially use the first order scheme for div(phi,nuTilda)
418 
419  p_, U_, phi_, etc: State variables in OpenFOAM
420 
421  Output:
422  nuTildaRes_: If solveTurbState_ == 0, update the residual field variable
423 
424  nuTilda_: If solveTurbState_ == 1, update nuTilda
425  */
426 
427  // Copy and modify based on the "correct" function
428 
429  word divNuTildaScheme = "div(phi,nuTilda)";
430 
431  label isPC = 0;
432 
433  label printToScreen = options.lookupOrDefault("printToScreen", 0);
434 
435  if (!solveTurbState_)
436  {
437  isPC = options.getLabel("isPC");
438 
439  if (isPC)
440  {
441  divNuTildaScheme = "div(pc)";
442  }
443  }
444 
445  const volScalarField chi(this->chi());
446  const volScalarField fv1(this->fv1(chi));
447 
448  const volScalarField Stilda(this->Stilda(chi, fv1));
449 
450  volScalarField rho = this->rho();
451 
452  tmp<fvScalarMatrix> nuTildaEqn(
453  fvm::ddt(phase_, rho, nuTilda_)
454  + fvm::div(phaseRhoPhi_, nuTilda_, divNuTildaScheme)
455  - fvm::laplacian(phase_ * rho * DnuTildaEff(), nuTilda_)
456  - Cb2_ / sigmaNut_ * phase_ * rho * magSqr(fvc::grad(nuTilda_))
458  - fvm::Sp(Cw1_ * phase_ * rho * fw(Stilda) * nuTilda_ / sqr(y_), nuTilda_));
459 
460  nuTildaEqn.ref().relax();
461 
462  if (solveTurbState_)
463  {
464  // get the solver performance info such as initial
465  // and final residuals
466  SolverPerformance<scalar> solverNuTilda = solve(nuTildaEqn);
467 
468  DAUtility::primalResidualControl(solverNuTilda, printToScreen, "nuTilda", daGlobalVar_.primalMaxRes);
469 
470  DAUtility::boundVar(allOptions_, nuTilda_, printToScreen);
471  nuTilda_.correctBoundaryConditions();
472 
473  // ***************** NOTE*****************
474  // In the original SA, it is correctNut(fv1) and fv1 is not
475  // updated based on the latest nuTilda. We use correctNut which
476  // recompute fv1 with the latest nuTilda
477  this->correctNut();
478  }
479  else
480  {
481  // calculate residuals
482  nuTildaRes_ = nuTildaEqn.ref() & nuTilda_;
483  // need to normalize residuals
484  normalizeResiduals(nuTildaRes);
485  }
486 
487  return;
488 }
489 
491  const word varName,
492  scalarField& diag,
493  scalarField& upper,
494  scalarField& lower)
495 {
496  /*
497  Description:
498  return the diag(), upper(), and lower() scalarFields from the turbulence model's fvMatrix
499  this will be use to compute the preconditioner matrix
500  */
501 
502  if (varName != "nuTilda")
503  {
504  FatalErrorIn(
505  "varName not valid. It has to be nuTilda")
506  << exit(FatalError);
507  }
508 
509  const volScalarField chi(this->chi());
510  const volScalarField fv1(this->fv1(chi));
511 
512  const volScalarField Stilda(this->Stilda(chi, fv1));
513 
514  volScalarField rho = this->rho();
515 
516  fvScalarMatrix nuTildaEqn(
517  fvm::ddt(phase_, rho, nuTilda_)
518  + fvm::div(phaseRhoPhi_, nuTilda_, "div(pc)")
519  - fvm::laplacian(phase_ * rho * DnuTildaEff(), nuTilda_)
520  - Cb2_ / sigmaNut_ * phase_ * rho * magSqr(fvc::grad(nuTilda_))
522  - fvm::Sp(Cw1_ * phase_ * rho * fw(Stilda) * nuTilda_ / sqr(y_), nuTilda_));
523 
524  nuTildaEqn.relax();
525 
526  diag = nuTildaEqn.D();
527  upper = nuTildaEqn.upper();
528  lower = nuTildaEqn.lower();
529 }
530 
532  const word varName,
533  const scalarList& rhs,
534  scalarList& dPsi)
535 {
536  if (varName != "nuTilda")
537  {
538  FatalErrorIn(
539  "varName not valid. It has to be nuTilda")
540  << exit(FatalError);
541  }
542 
543  const fvSolution& myFvSolution = mesh_.thisDb().lookupObject<fvSolution>("fvSolution");
544  dictionary solverDictNuTilda = myFvSolution.subDict("solvers").subDict("nuTilda");
545 
546  // solve the fvMatrixT field with given rhs and solution
547  if (psiNuTildaPC_.empty())
548  {
550  const volScalarField chi(this->chi());
551  const volScalarField fv1(this->fv1(chi));
552 
553  const volScalarField Stilda(this->Stilda(chi, fv1));
554 
555  volScalarField rho = this->rho();
556 
557  // for the PC, we need only fvm:: terms
558  psiNuTildaPC_.reset(new fvScalarMatrix(
559  fvm::ddt(phase_, rho, dPsiNuTilda_)
560  + fvm::div(phaseRhoPhi_, dPsiNuTilda_, "div(phi,nuTilda)")
561  - fvm::laplacian(phase_ * rho * DnuTildaEff(), dPsiNuTilda_)
562  - Cb2_ / sigmaNut_ * phase_ * rho * magSqr(fvc::grad(dPsiNuTilda_))
563  + fvm::Sp(Cw1_ * phase_ * rho * fw(Stilda) * dPsiNuTilda_ / sqr(y_), dPsiNuTilda_)));
564 
565  psiNuTildaPC_->relax();
566 
567  DAUtility::swapLists<scalar>(psiNuTildaPC_->upper(), psiNuTildaPC_->lower());
568 
569  // make sure the boundary contribution to source hrs is zero
570  forAll(dPsiNuTilda_.boundaryField(), patchI)
571  {
572  forAll(dPsiNuTilda_.boundaryField()[patchI], faceI)
573  {
574  psiNuTildaPC_->boundaryCoeffs()[patchI][faceI] = 0.0;
575  }
576  }
577  }
578 
579  // force to use zeros as the initial guess
580  forAll(dPsiNuTilda_, cellI)
581  {
582  dPsiNuTilda_[cellI] = 0.0;
583  }
584 
585  // set the rhs
586  forAll(psiNuTildaPC_->source(), cellI)
587  {
588  psiNuTildaPC_->source()[cellI] = rhs[cellI];
589  }
590 
591  // solve
592  psiNuTildaPC_->solve(solverDictNuTilda);
593 
594  // return the solution
595  forAll(dPsiNuTilda_, cellI)
596  {
597  dPsi[cellI] = dPsiNuTilda_[cellI];
598  }
599 }
600 
601 void DASpalartAllmaras::getTurbProdOverDestruct(volScalarField& PoD) const
602 {
603  /*
604  Description:
605  Return the value of the production over destruction term from the turbulence model
606  */
607 
608  const volScalarField chi(this->chi());
609  const volScalarField fv1(this->fv1(chi));
610 
611  const volScalarField Stilda(this->Stilda(chi, fv1));
612 
613  volScalarField rho = this->rho();
614 
615  volScalarField P = Cb1_ * phase_ * rho * Stilda * nuTilda_;
616  volScalarField D = Cw1_ * phase_ * rho * fw(Stilda) * sqr(nuTilda_ / y_);
617 
618  forAll(P, cellI)
619  {
620  PoD[cellI] = P[cellI] / (D[cellI] + P[cellI] + 1e-16);
621  }
622 }
623 
624 void DASpalartAllmaras::getTurbConvOverProd(volScalarField& CoP) const
625 {
626  /*
627  Description:
628  Return the value of the convective over production term from the turbulence model
629  */
630 
631  const volScalarField chi(this->chi());
632  const volScalarField fv1(this->fv1(chi));
633 
634  const volScalarField Stilda(this->Stilda(chi, fv1));
635 
636  volScalarField rho = this->rho();
637 
638  volScalarField P = Cb1_ * phase_ * rho * Stilda * nuTilda_;
639  volScalarField C = fvc::div(phaseRhoPhi_, nuTilda_);
640 
641  forAll(P, cellI)
642  {
643  CoP[cellI] = C[cellI] / (P[cellI] + C[cellI] + 1e-16);
644  }
645 }
646 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
647 
648 } // End namespace Foam
649 
650 // ************************************************************************* //
Foam::DASpalartAllmaras::correct
virtual void correct(label printToScreen)
solve the residual equations and update the state
Definition: DASpalartAllmaras.C:386
Foam::DASpalartAllmaras::sigmaNut_
dimensionedScalar sigmaNut_
Definition: DASpalartAllmaras.H:55
solve
nuTilda1Eqn solve(solverDictNuTilda)
D
volVectorField & D
Definition: createRefsSolidDisplacement.H:8
DASpalartAllmaras.H
Foam::DASpalartAllmaras::Cw1_
dimensionedScalar Cw1_
Definition: DASpalartAllmaras.H:59
Foam::DASpalartAllmaras::DASpalartAllmaras
DASpalartAllmaras(const word modelType, const fvMesh &mesh, const DAOption &daOption)
Definition: DASpalartAllmaras.C:41
Foam::DASpalartAllmaras::solveTurbState_
label solveTurbState_
whether to solve for turb states
Definition: DASpalartAllmaras.H:97
Foam::DASpalartAllmaras::Cv1_
dimensionedScalar Cv1_
Definition: DASpalartAllmaras.H:62
Foam::DAOption
Definition: DAOption.H:29
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFunction, DAFunctionForce, dictionary)
Foam::DASpalartAllmaras::correctStateResidualModelCon
virtual void correctStateResidualModelCon(List< List< word >> &stateCon) const
update the original variable connectivity for the adjoint state residuals in stateCon
Definition: DASpalartAllmaras.C:257
Foam::DASpalartAllmaras::getTurbProdOverDestruct
virtual void getTurbProdOverDestruct(volScalarField &PoD) const
return the value of the destruction term from the turbulence model
Definition: DASpalartAllmaras.C:601
normalizeResiduals
#define normalizeResiduals(resName)
Definition: DAMacroFunctions.H:28
Foam::DASpalartAllmaras::betaFINuTilda_
volScalarField betaFINuTilda_
beta field for field inversion
Definition: DASpalartAllmaras.H:94
Foam::DASpalartAllmaras::Cb1_
dimensionedScalar Cb1_
Definition: DASpalartAllmaras.H:57
Foam::DATurbulenceModel::phaseRhoPhi_
surfaceScalarField & phaseRhoPhi_
phase*phi*density field
Definition: DATurbulenceModel.H:89
Foam::DATurbulenceModel::allOptions_
const dictionary & allOptions_
all DAFoam options
Definition: DATurbulenceModel.H:71
Foam::DATurbulenceModel::rho
tmp< volScalarField > rho() const
get the density field
Definition: DATurbulenceModel.C:294
Foam::DASpalartAllmaras::Cw2_
dimensionedScalar Cw2_
Definition: DASpalartAllmaras.H:60
Foam::DASpalartAllmaras::kappa_
dimensionedScalar kappa_
Definition: DASpalartAllmaras.H:56
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
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:475
Foam::DAUtility::primalResidualControl
static void primalResidualControl(const SolverPerformance< scalar > &solverP, const label printToScreen, const word varName, scalar &primalMaxRes)
control when to print the residual and also compute the maxInitRes
Definition: DAUtility.C:734
Foam::DATurbulenceModel::daGlobalVar_
DAGlobalVar & daGlobalVar_
global variable
Definition: DATurbulenceModel.H:74
Foam::DATurbulenceModel::U_
volVectorField & U_
velocity
Definition: DATurbulenceModel.H:80
Foam::DATurbulenceModel::nut_
volScalarField & nut_
turbulence viscosity
Definition: DATurbulenceModel.H:77
Foam::DAGlobalVar::primalMaxRes
scalar primalMaxRes
the maximal residual for primal this var will be used in multiple classes
Definition: DAGlobalVar.H:62
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFunction, 0)
Foam::DASpalartAllmaras::Cs_
dimensionedScalar Cs_
Definition: DASpalartAllmaras.H:63
Foam::DASpalartAllmaras::psiNuTildaPC_
autoPtr< fvScalarMatrix > psiNuTildaPC_
Definition: DASpalartAllmaras.H:99
Foam
Definition: checkGeometry.C:32
Foam::DASpalartAllmaras::chi
tmp< volScalarField > chi() const
Definition: DASpalartAllmaras.C:124
Foam::DASpalartAllmaras::getFvMatrixFields
virtual void getFvMatrixFields(const word varName, scalarField &diag, scalarField &upper, scalarField &lower)
return the diag(), upper(), and lower() scalarFields from the turbulence model's fvMatrix
Definition: DASpalartAllmaras.C:490
Foam::DASpalartAllmaras::calcResiduals
virtual void calcResiduals(const dictionary &options)
compute the turbulence residuals
Definition: DASpalartAllmaras.C:407
Foam::DASpalartAllmaras::Stilda
tmp< volScalarField > Stilda(const volScalarField &chi, const volScalarField &fv1) const
Definition: DASpalartAllmaras.C:143
Foam::DASpalartAllmaras::dPsiNuTilda_
volScalarField dPsiNuTilda_
Definition: DASpalartAllmaras.H:101
Foam::DATurbulenceModel::nu
tmp< volScalarField > nu() const
get the nu field
Definition: DATurbulenceModel.C:267
forAll
forAll(nuTilda1, cellI)
Definition: nuTilda1EqnIrkPimple.H:19
Foam::DASpalartAllmaras::nuTildaRes_
volScalarField nuTildaRes_
Definition: DASpalartAllmaras.H:87
Foam::DASpalartAllmaras::addModelResidualCon
virtual void addModelResidualCon(HashTable< List< List< word >>> &allCon) const
add the model residual connectivity to stateCon
Definition: DASpalartAllmaras.C:307
Foam::DATurbulenceModel::correctAlphat
void correctAlphat()
update alphat
Definition: DATurbulenceModel.C:203
Foam::DASpalartAllmaras::nuTilda_
volScalarField & nuTilda_
Definition: DASpalartAllmaras.H:86
Foam::DATurbulenceModel
Definition: DATurbulenceModel.H:52
Foam::DASpalartAllmaras::Cw3_
dimensionedScalar Cw3_
Definition: DASpalartAllmaras.H:61
Foam::DASpalartAllmaras::correctBoundaryConditions
virtual void correctBoundaryConditions()
update turbulence variable boundary values
Definition: DASpalartAllmaras.C:235
Foam::DATurbulenceModel::mesh_
const fvMesh & mesh_
fvMesh
Definition: DATurbulenceModel.H:65
Foam::DASpalartAllmaras::correctNut
virtual void correctNut()
update nut based on other turbulence variables and update the BCs
Definition: DASpalartAllmaras.C:215
Foam::DASpalartAllmaras::fv1
tmp< volScalarField > fv1(const volScalarField &chi) const
Definition: DASpalartAllmaras.C:129
Foam::DASpalartAllmaras::DnuTildaEff
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
Definition: DASpalartAllmaras.C:174
Foam::DASpalartAllmaras::correctModelStates
virtual void correctModelStates(wordList &modelStates) const
update the turbulence state for DAStateInfo::regStates_
Definition: DASpalartAllmaras.C:181
Foam::DASpalartAllmaras::solveAdjointFP
virtual void solveAdjointFP(const word varName, const scalarList &rhs, scalarList &dPsi)
solve the fvMatrixT field with given rhs and solution
Definition: DASpalartAllmaras.C:531
Foam::DASpalartAllmaras::Cb2_
dimensionedScalar Cb2_
Definition: DASpalartAllmaras.H:58
Foam::DASpalartAllmaras::updateIntermediateVariables
virtual void updateIntermediateVariables()
update any intermediate variables that are dependent on state variables and are used in calcResiduals
Definition: DASpalartAllmaras.C:246
Foam::DASpalartAllmaras::fv2
tmp< volScalarField > fv2(const volScalarField &chi, const volScalarField &fv1) const
Definition: DASpalartAllmaras.C:136
Foam::DASpalartAllmaras::fw
tmp< volScalarField > fw(const volScalarField &Stilda) const
Definition: DASpalartAllmaras.C:156
Foam::DASpalartAllmaras::getTurbConvOverProd
virtual void getTurbConvOverProd(volScalarField &CoP) const
return the value of the convective over production term from the turbulence model
Definition: DASpalartAllmaras.C:624
Foam::DATurbulenceModel::turbModelType_
word turbModelType_
whether the turbulence model is incompressible or compressible
Definition: DATurbulenceModel.H:95
Foam::DATurbulenceModel::phase_
volScalarField phase_
phase field
Definition: DATurbulenceModel.H:86
Foam::DASpalartAllmaras::y_
const volScalarField & y_
3D wall distance
Definition: DASpalartAllmaras.H:91