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