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