DASpalartAllmarasFv3.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/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 "DASpalartAllmarasFv3.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 defineTypeNameAndDebug(DASpalartAllmarasFv3, 0);
38 addToRunTimeSelectionTable(DATurbulenceModel, DASpalartAllmarasFv3, 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  Cv2_(dimensioned<scalar>::lookupOrAddToDict(
78  "Cv2",
79  this->coeffDict_,
80  5.0)),
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 #ifdef CompressibleFlow
93  dimensionedScalar("nuTildaRes", dimensionSet(1, -1, -2, 0, 0, 0, 0), 0.0),
94 #endif
95 #ifdef IncompressibleFlow
96  dimensionedScalar("nuTildaRes", dimensionSet(0, 2, -2, 0, 0, 0, 0), 0.0),
97 #endif
98  zeroGradientFvPatchField<scalar>::typeName),
99  // pseudoNuTilda_ and pseudoNuTildaEqn_ for solving adjoint equation
100  pseudoNuTilda_(
101  IOobject(
102  "pseudoNuTilda",
103  mesh.time().timeName(),
104  mesh,
105  IOobject::NO_READ,
106  IOobject::NO_WRITE),
107  nuTilda_),
108  pseudoNuTildaEqn_(fvm::div(phi_, pseudoNuTilda_, "div(phi,nuTilda)")),
109  y_(mesh.thisDb().lookupObject<volScalarField>("yWall"))
110 {
111 
112  // initialize printInterval_ we need to check whether it is a steady state
113  // or unsteady primal solver
114  IOdictionary fvSchemes(
115  IOobject(
116  "fvSchemes",
117  mesh.time().system(),
118  mesh,
119  IOobject::MUST_READ,
120  IOobject::NO_WRITE,
121  false));
122  word ddtScheme = word(fvSchemes.subDict("ddtSchemes").lookup("default"));
123  if (ddtScheme == "steadyState")
124  {
126  daOption.getAllOptions().lookupOrDefault<label>("printInterval", 100);
127  }
128  else
129  {
131  daOption.getAllOptions().lookupOrDefault<label>("printIntervalUnsteady", 500);
132  }
133 
134  // get fvSolution and fvSchemes info for fixed-point adjoint
135  const fvSolution& myFvSolution = mesh.thisDb().lookupObject<fvSolution>("fvSolution");
136  solverDictNuTilda_ = myFvSolution.subDict("solvers").subDict("nuTilda");
137  if (myFvSolution.found("relaxationFactors"))
138  {
139  if (myFvSolution.subDict("relaxationFactors").found("equations"))
140  {
141  if (myFvSolution.subDict("relaxationFactors").subDict("equations").found("nuTilda"))
142  {
143  relaxNuTildaEqn_ = myFvSolution.subDict("relaxationFactors").subDict("equations").getScalar("nuTilda");
144  }
145  }
146  }
147 }
148 
149 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
150 
151 // SA member functions. these functions are copied from
152 tmp<volScalarField> DASpalartAllmarasFv3::chi() const
153 {
154  return nuTilda_ / this->nu();
155 }
156 
157 tmp<volScalarField> DASpalartAllmarasFv3::fv1(
158  const volScalarField& chi) const
159 {
160  const volScalarField chi3(pow3(chi));
161  return chi3 / (chi3 + pow3(Cv1_));
162 }
163 
164 tmp<volScalarField> DASpalartAllmarasFv3::fv2(
165  const volScalarField& chi,
166  const volScalarField& fv1) const
167 {
168  return 1.0 / pow3(scalar(1) + chi / Cv2_);
169 }
170 
171 tmp<volScalarField> DASpalartAllmarasFv3::fv3(
172  const volScalarField& chi,
173  const volScalarField& fv1) const
174 {
175 
176  const volScalarField chiByCv2((1 / Cv2_) * chi);
177 
178  return (scalar(1) + chi * fv1)
179  * (1 / Cv2_)
180  * (3 * (scalar(1) + chiByCv2) + sqr(chiByCv2))
181  / pow3(scalar(1) + chiByCv2);
182 }
183 
184 tmp<volScalarField> DASpalartAllmarasFv3::fw(
185  const volScalarField& Stilda) const
186 {
187  volScalarField r(
188  min(
189  nuTilda_
190  / (max(
191  Stilda,
192  dimensionedScalar("SMALL", Stilda.dimensions(), SMALL))
193  * sqr(kappa_ * y_)),
194  scalar(10.0)));
195  r.boundaryFieldRef() == 0.0;
196 
197  const volScalarField g(r + Cw2_ * (pow6(r) - r));
198 
199  return g * pow((1.0 + pow6(Cw3_)) / (pow6(g) + pow6(Cw3_)), 1.0 / 6.0);
200 }
201 
202 tmp<volScalarField> DASpalartAllmarasFv3::DnuTildaEff() const
203 {
204  return tmp<volScalarField>(
205  new volScalarField("DnuTildaEff", (nuTilda_ + this->nu()) / sigmaNut_));
206 }
207 
208 // Augmented functions
209 void DASpalartAllmarasFv3::correctModelStates(wordList& modelStates) const
210 {
211  /*
212  Description:
213  Update the name in modelStates based on the selected physical model at runtime
214 
215  Example:
216  In DAStateInfo, if the modelStates reads:
217 
218  modelStates = {"nut"}
219 
220  then for the SA model, calling correctModelStates(modelStates) will give:
221 
222  modelStates={"nuTilda"}
223 
224  while calling correctModelStates(modelStates) for the SST model will give
225 
226  modelStates={"k","omega"}
227 
228  We don't udpate the names for the radiation model becasue users are
229  supposed to set modelStates={"G"}
230  */
231 
232  // replace nut with nuTilda
233  forAll(modelStates, idxI)
234  {
235  word stateName = modelStates[idxI];
236  if (stateName == "nut")
237  {
238  modelStates[idxI] = "nuTilda";
239  }
240  }
241 }
242 
244 {
245  /*
246  Description:
247  Update nut based on other turbulence variables and update the BCs
248  Also update alphat if is present
249  */
250 
251  const volScalarField chi(this->chi());
252  const volScalarField fv1(this->fv1(chi));
253  nut_ = nuTilda_ * fv1;
254 
255  nut_.correctBoundaryConditions();
256 
257  // this is basically BasicTurbulenceModel::correctNut();
258  this->correctAlphat();
259 
260  return;
261 }
262 
264 {
265  /*
266  Description:
267  Update turbulence variable boundary values
268  */
269 
270  // correct the BCs for the perturbed fields
271  nuTilda_.correctBoundaryConditions();
272 }
273 
275 {
276  /*
277  Description:
278  Update nut based on nuTilda. Note: we need to update nut and its BC since we
279  may have perturbed other turbulence vars that affect the nut values
280  */
281 
282  this->correctNut();
283 }
284 
285 void DASpalartAllmarasFv3::correctStateResidualModelCon(List<List<word>>& stateCon) const
286 {
287  /*
288  Description:
289  Update the original variable connectivity for the adjoint state
290  residuals in stateCon. Basically, we modify/add state variables based on the
291  original model variables defined in stateCon.
292 
293  Input:
294 
295  stateResCon: the connectivity levels for a state residual, defined in Foam::DAJacCon
296 
297  Example:
298  If stateCon reads:
299  stateCon=
300  {
301  {"U", "p", "nut"},
302  {"p"}
303  }
304 
305  For the SA turbulence model, calling this function for will get a new stateCon
306  stateCon=
307  {
308  {"U", "p", "nuTilda"},
309  {"p"}
310  }
311 
312  For the SST turbulence model, calling this function will give
313  stateCon=
314  {
315  {"U", "p", "k", "omega"},
316  {"p", "U"}
317  }
318  ***NOTE***: we add a extra level of U connectivity because nut is
319  related to grad(U), k, and omega in SST!
320  */
321 
322  forAll(stateCon, idxI)
323  {
324  forAll(stateCon[idxI], idxJ)
325  {
326  word conStateName = stateCon[idxI][idxJ];
327  if (conStateName == "nut")
328  {
329  stateCon[idxI][idxJ] = "nuTilda";
330  }
331  }
332  }
333 }
334 
335 void DASpalartAllmarasFv3::addModelResidualCon(HashTable<List<List<word>>>& allCon) const
336 {
337  /*
338  Description:
339  Add the connectivity levels for all physical model residuals to allCon
340 
341  Input:
342  allCon: the connectivity levels for all state residual, defined in DAJacCon
343 
344  Example:
345  If stateCon reads:
346  allCon=
347  {
348  "URes":
349  {
350  {"U", "p", "nut"},
351  {"p"}
352  }
353  }
354 
355  For the SA turbulence model, calling this function for will get a new stateCon,
356  something like this:
357  allCon=
358  {
359  "URes":
360  {
361  {"U", "p", "nuTilda"},
362  {"p"}
363  },
364  "nuTildaRes":
365  {
366  {"U", "phi", "nuTilda"},
367  {"U"}
368  }
369  }
370 
371  */
372 
373  word pName;
374 
375  if (mesh_.thisDb().foundObject<volScalarField>("p"))
376  {
377  pName = "p";
378  }
379  else if (mesh_.thisDb().foundObject<volScalarField>("p_rgh"))
380  {
381  pName = "p_rgh";
382  }
383  else
384  {
385  FatalErrorIn(
386  "Neither p nor p_rgh was found in mesh.thisDb()!"
387  "addModelResidualCon failed to setup turbulence residuals!")
388  << exit(FatalError);
389  }
390 
391  // NOTE: for compressible flow, it depends on rho so we need to add T and p
392 #ifdef IncompressibleFlow
393  allCon.set(
394  "nuTildaRes",
395  {
396  {"U", "nuTilda", "phi"}, // lv0
397  {"U", "nuTilda"}, // lv1
398  {"nuTilda"} // lv2
399  });
400 #endif
401 
402 #ifdef CompressibleFlow
403  allCon.set(
404  "nuTildaRes",
405  {
406  {"U", "T", pName, "nuTilda", "phi"}, // lv0
407  {"U", "T", pName, "nuTilda"}, // lv1
408  {"T", pName, "nuTilda"} // lv2
409  });
410 #endif
411 }
412 
414 {
415  /*
416  Descroption:
417  Solve the residual equations and update the state. This function will be called
418  by the DASolver. It is needed because we want to control the output frequency
419  of the residual convergence every 100 steps. If using the correct from turbulence
420  it will output residual every step which will be too much of information.
421  */
422 
423  // We set the flag solveTurbState_ to 1 such that in the calcResiduals function
424  // we will solve and update nuTilda
425  solveTurbState_ = 1;
426  dictionary dummyOptions;
427  this->calcResiduals(dummyOptions);
428  // after it, we reset solveTurbState_ = 0 such that calcResiduals will not
429  // update nuTilda when calling from the adjoint class, i.e., solveAdjoint from DASolver.
430  solveTurbState_ = 0;
431 }
432 
433 void DASpalartAllmarasFv3::calcResiduals(const dictionary& options)
434 {
435  /*
436  Descroption:
437  If solveTurbState_ == 1, this function solve and update nuTilda, and
438  is the same as calling turbulence.correct(). If solveTurbState_ == 0,
439  this function compute residuals for turbulence variables, e.g., nuTildaRes_
440 
441  Input:
442  options.isPC: 1 means computing residuals for preconditioner matrix.
443  This essentially use the first order scheme for div(phi,nuTilda)
444 
445  p_, U_, phi_, etc: State variables in OpenFOAM
446 
447  Output:
448  nuTildaRes_: If solveTurbState_ == 0, update the residual field variable
449 
450  nuTilda_: If solveTurbState_ == 1, update nuTilda
451  */
452 
453  // Copy and modify based on the "correct" function
454 
455  label printToScreen = this->isPrintTime(mesh_.time(), printInterval_);
456 
457  word divNuTildaScheme = "div(phi,nuTilda)";
458 
459  label isPC = 0;
460 
461  if (!solveTurbState_)
462  {
463  isPC = options.getLabel("isPC");
464 
465  if (isPC)
466  {
467  divNuTildaScheme = "div(pc)";
468  }
469  }
470 
471  const volScalarField chi(this->chi());
472  const volScalarField fv1(this->fv1(chi));
473 
474  const volScalarField Stilda(
475  this->fv3(chi, fv1) * ::sqrt(2.0) * mag(skew(fvc::grad(U_)))
476  + this->fv2(chi, fv1) * nuTilda_ / sqr(kappa_ * y_));
477 
478  tmp<fvScalarMatrix> nuTildaEqn(
479  fvm::ddt(phase_, rho_, nuTilda_)
480  + fvm::div(phaseRhoPhi_, nuTilda_, divNuTildaScheme)
481  - fvm::laplacian(phase_ * rho_ * DnuTildaEff(), nuTilda_)
482  - Cb2_ / sigmaNut_ * phase_ * rho_ * magSqr(fvc::grad(nuTilda_))
483  == Cb1_ * phase_ * rho_ * Stilda * nuTilda_
484  - fvm::Sp(Cw1_ * phase_ * rho_ * fw(Stilda) * nuTilda_ / sqr(y_), nuTilda_));
485 
486  nuTildaEqn.ref().relax();
487 
488  if (solveTurbState_)
489  {
490 
491  // get the solver performance info such as initial
492  // and final residuals
493  SolverPerformance<scalar> solverNuTilda = solve(nuTildaEqn);
494  if (printToScreen)
495  {
496  Info << "nuTilda Initial residual: " << solverNuTilda.initialResidual() << endl
497  << " Final residual: " << solverNuTilda.finalResidual() << endl;
498  }
499 
500  DAUtility::boundVar(allOptions_, nuTilda_, printToScreen);
501  nuTilda_.correctBoundaryConditions();
502 
503  // ***************** NOTE*****************
504  // In the original SA, it is correctNut(fv1) and fv1 is not
505  // updated based on the latest nuTilda. We use correctNut which
506  // recompute fv1 with the latest nuTilda
507  this->correctNut();
508  }
509  else
510  {
511  // calculate residuals
512  nuTildaRes_ = nuTildaEqn.ref() & nuTilda_;
513  // need to normalize residuals
514  normalizeResiduals(nuTildaRes);
515  }
516 
517  return;
518 }
519 
521  const volScalarField& mySource,
522  volScalarField& pseudoNuTilda)
523 {
524  /*
525  Description:
526  Inverse transpose product, M_nuTilda^(-T)
527  Based on inverseProduct_nuTildaEqn from simpleFoamSAPrimal, but swaping upper() and lower()
528  We won't ADR this function, so we can treat most of the arguments as const
529  */
530 
531  // Make sure pseudoNuTilda = nuTilda;
532  //if (pseudoNuTildaEqnInitialized_ == 0)
533  //{
534  /*
535  pseudoNuTilda_ = nuTilda_;
536 
537  const volScalarField chi(this->chi());
538  const volScalarField fv1(this->fv1(chi));
539 
540  // Get myStilda
541  const volScalarField Stilda(
542  this->fv3(chi, fv1) * ::sqrt(2.0) * mag(skew(fvc::grad(U_)))
543  + this->fv2(chi, fv1) * nuTilda_ / sqr(kappa_ * y_));
544 
545  // Get the pseudoNuTildaEqn,
546  // the most important thing here is to make sure the l.h.s. mathces that of nuTildaEqn.
547  // Some explicit terms that only contributes to the r.h.s. are diabled
548  pseudoNuTildaEqn_ =
549  //fvm::ddt(pseudoNuTilda_)
550  fvm::div(phi_, pseudoNuTilda_, "div(phi,nuTilda)")
551  - fvm::laplacian(DnuTildaEff(), pseudoNuTilda_)
552  + fvm::Sp(Cw1_ * fw(Stilda) * pseudoNuTilda_ / sqr(y_), pseudoNuTilda_);
553  pseudoNuTildaEqn_.relax(relaxNuTildaEqn_);
554 
555  // Swap upper() and lower()
556  List<scalar> temp = pseudoNuTildaEqn_.upper();
557  pseudoNuTildaEqn_.upper() = pseudoNuTildaEqn_.lower();
558  pseudoNuTildaEqn_.lower() = temp;
559 
560  // mark it as initialized
561  // pseudoNuTildaEqnInitialized_ = 1;
562  //}
563 
564  // Overwrite the r.h.s.
565  pseudoNuTildaEqn_.source() = mySource.primitiveField();
566 
567  // Make sure that boundary contribution to source is zero,
568  // Alternatively, we can deduct source by boundary contribution, so that it would cancel out during solve.
569  forAll(pseudoNuTilda_.boundaryField(), patchI)
570  {
571  const fvPatch& pp = pseudoNuTilda_.boundaryField()[patchI].patch();
572  forAll(pp, faceI)
573  {
574  label cellI = pp.faceCells()[faceI];
575  pseudoNuTildaEqn_.source()[cellI] -= pseudoNuTildaEqn_.boundaryCoeffs()[patchI][faceI];
576  }
577  }
578 
579  // Before solve, force xxEqn.psi to be solved into all zero
580  // This ensures the zero (internal) initial guess
581  forAll(pseudoNuTilda_.primitiveFieldRef(), cellI)
582  {
583  pseudoNuTilda_.primitiveFieldRef()[cellI] = 0;
584  }
585  // Solve using the zero (internal) initial guess
586  pseudoNuTildaEqn_.solve(solverDictNuTilda_);
587 
588  forAll(pseudoNuTilda, cellI)
589  {
590  pseudoNuTilda[cellI] = pseudoNuTilda_[cellI];
591  }
592 */
593 }
594 
596 {
597  /*
598  Description:
599  construct the pseudo nuTildaEqn,
600  which is nuTildaEqn with the lhs upper and lower arrays swapped,
601  we also don't care about the rhs of pseudo nuTildaEqn.
602  */
603 
604  // Make sure pseudoNuTilda is indeed identical to nuTilda
606  pseudoNuTilda_.correctBoundaryConditions();
607 
608  const volScalarField chi(this->chi());
609  const volScalarField fv1(this->fv1(chi));
610 
611  // Get myStilda
612  const volScalarField Stilda(
613  this->fv3(chi, fv1) * ::sqrt(2.0) * mag(skew(fvc::grad(U_)))
614  + this->fv2(chi, fv1) * nuTilda_ / sqr(kappa_ * y_));
615 
616  // Get the pseudoNuTildaEqn,
617  // the most important thing here is to make sure the l.h.s. mathces that of nuTildaEqn.
618  // Some explicit terms that only contributes to the r.h.s. are diabled
620  //fvm::ddt(pseudoNuTilda_)
621  fvm::div(phi_, pseudoNuTilda_, "div(phi,nuTilda)")
622  - fvm::laplacian(DnuTildaEff(), pseudoNuTilda_)
623  + fvm::Sp(Cw1_ * fw(Stilda) * pseudoNuTilda_ / sqr(y_), pseudoNuTilda_);
625 
626  // Swap upper() and lower()
627  // We will use the swap function once it's being moved to the DAUtility
628  List<scalar> temp = pseudoNuTildaEqn_.upper();
629  pseudoNuTildaEqn_.upper() = pseudoNuTildaEqn_.lower();
630  pseudoNuTildaEqn_.lower() = temp;
631 }
632 
633 void DASpalartAllmarasFv3::rhsSolvePseudoNuTildaEqn(const volScalarField& nuTildaSource)
634 {
635  /*
636  Description:
637  solve the pseudo nuTildaEqn with a overwritten rhs
638  */
639 
640  // Overwrite the r.h.s.
641  pseudoNuTildaEqn_.source() = nuTildaSource.primitiveField();
642 
643  // Make sure that boundary contribution to source is zero,
644  // Alternatively, we can deduct source by boundary contribution, so that it would cancel out during solve.
645  forAll(pseudoNuTilda_.boundaryField(), patchI)
646  {
647  const fvPatch& pp = pseudoNuTilda_.boundaryField()[patchI].patch();
648  forAll(pp, faceI)
649  {
650  label cellI = pp.faceCells()[faceI];
651  pseudoNuTildaEqn_.source()[cellI] -= pseudoNuTildaEqn_.boundaryCoeffs()[patchI][faceI];
652  }
653  }
654 
655  // Before solve, force xxEqn.psi to be solved into all zero
656  // This ensures the zero (internal) initial guess
657  forAll(pseudoNuTilda_.primitiveFieldRef(), cellI)
658  {
659  pseudoNuTilda_.primitiveFieldRef()[cellI] = 0;
660  }
661  // Solve using the zero (internal) initial guess
663 }
664 
665 void DASpalartAllmarasFv3::calcLduResidualTurb(volScalarField& nuTildaRes)
666 {
667  /*
668  Description:
669  calculate the turbulence residual using LDU matrix
670  */
671 
672  const volScalarField chi(this->chi());
673  const volScalarField fv1(this->fv1(chi));
674 
675  // Get myStilda
676  const volScalarField Stilda(
677  this->fv3(chi, fv1) * ::sqrt(2.0) * mag(skew(fvc::grad(U_)))
678  + this->fv2(chi, fv1) * nuTilda_ / sqr(kappa_ * y_));
679 
680  // Construct nuTildaEqn using our own SA implementation
681  fvScalarMatrix nuTildaEqn(
682  fvm::div(phi_, nuTilda_, "div(phi,nuTilda)")
683  - fvm::laplacian(DnuTildaEff(), nuTilda_)
684  - Cb2_ / sigmaNut_ * magSqr(fvc::grad(nuTilda_))
685  == Cb1_ * Stilda * nuTilda_
686  - fvm::Sp(Cw1_ * fw(Stilda) * nuTilda_ / sqr(y_), nuTilda_));
687 
688  List<scalar>& nuTildaSource = nuTildaEqn.source();
689  List<scalar>& nuTildaDiag = nuTildaEqn.diag();
690 
691  // Initiate nuTildaRes, with no boundary contribution
692  for (label i = 0; i < nuTilda_.size(); i++)
693  {
694  nuTildaRes[i] = nuTildaDiag[i] * nuTilda_[i] - nuTildaSource[i];
695  }
696  nuTildaRes.primitiveFieldRef() -= nuTildaEqn.lduMatrix::H(nuTilda_);
697 
698  // Boundary correction
699  forAll(nuTilda_.boundaryField(), patchI)
700  {
701  const fvPatch& pp = nuTilda_.boundaryField()[patchI].patch();
702  forAll(pp, faceI)
703  {
704  // Both ways of getting cellI work
705  // Below is the previous way of getting the address
706  label cellI = pp.faceCells()[faceI];
707  // Below is using lduAddr().patchAddr(patchi)
708  //label cellI = nuTildaEqn.lduAddr().patchAddr(patchI)[faceI];
709  //myDiag[cellI] += TEqn.internalCoeffs()[patchI][faceI];
710  nuTildaRes[cellI] += nuTildaEqn.internalCoeffs()[patchI][faceI] * nuTilda_[cellI];
711  nuTildaRes[cellI] -= nuTildaEqn.boundaryCoeffs()[patchI][faceI];
712  }
713  }
714 
715  // Below is not necessary, but it doesn't hurt
716  nuTildaRes.correctBoundaryConditions();
717 }
718 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
719 
720 } // End namespace Foam
721 
722 // ************************************************************************* //
Foam::DASpalartAllmarasFv3::solverDictNuTilda_
dictionary solverDictNuTilda_
Definition: DASpalartAllmarasFv3.H:105
Foam::DASpalartAllmarasFv3::y_
const volScalarField & y_
3D wall distance
Definition: DASpalartAllmarasFv3.H:95
Foam::DASpalartAllmarasFv3::fv1
tmp< volScalarField > fv1(const volScalarField &chi) const
Definition: DASpalartAllmarasFv3.C:157
forAll
forAll(pseudoP.boundaryField(), patchI)
Definition: solvePseudoPEqn.H:10
Foam::DATurbulenceModel::phi_
surfaceScalarField & phi_
face flux
Definition: DATurbulenceModel.H:89
Foam::DASpalartAllmarasFv3::calcResiduals
virtual void calcResiduals(const dictionary &options)
compute the turbulence residuals
Definition: DASpalartAllmarasFv3.C:433
Foam::DAOption
Definition: DAOption.H:29
Foam::DASpalartAllmarasFv3::nuTilda_
volScalarField & nuTilda_
Definition: DASpalartAllmarasFv3.H:86
Foam::DASpalartAllmarasFv3::Cw3_
dimensionedScalar Cw3_
Definition: DASpalartAllmarasFv3.H:61
daOption
DAOption daOption(mesh, pyOptions_)
pseudoNuTilda
volScalarField & pseudoNuTilda
Definition: pseudoEqns.H:34
Foam::DASpalartAllmarasFv3::solveTurbState_
label solveTurbState_
whether to solve for turb states
Definition: DASpalartAllmarasFv3.H:98
normalizeResiduals
#define normalizeResiduals(resName)
Definition: DAMacroFunctions.H:28
DASpalartAllmarasFv3.H
Foam::DASpalartAllmarasFv3::Cw1_
dimensionedScalar Cw1_
Definition: DASpalartAllmarasFv3.H:59
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::DASpalartAllmarasFv3::correctModelStates
virtual void correctModelStates(wordList &modelStates) const
update the turbulence state for DAStateInfo::regStates_
Definition: DASpalartAllmarasFv3.C:209
Foam::DASpalartAllmarasFv3::relaxNuTildaEqn_
scalar relaxNuTildaEqn_
Definition: DASpalartAllmarasFv3.H:104
Foam::DASpalartAllmarasFv3::nuTildaRes_
volScalarField nuTildaRes_
Definition: DASpalartAllmarasFv3.H:87
Foam::DASpalartAllmarasFv3::Cb1_
dimensionedScalar Cb1_
Definition: DASpalartAllmarasFv3.H:57
Foam::DASpalartAllmarasFv3::Cv1_
dimensionedScalar Cv1_
Definition: DASpalartAllmarasFv3.H:62
Foam::DASpalartAllmarasFv3::correctNut
virtual void correctNut()
update nut based on other turbulence variables and update the BCs
Definition: DASpalartAllmarasFv3.C:243
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DASpalartAllmarasFv3::fv3
tmp< volScalarField > fv3(const volScalarField &chi, const volScalarField &fv1) const
Definition: DASpalartAllmarasFv3.C:171
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::DASpalartAllmarasFv3::pseudoNuTildaEqn_
fvScalarMatrix pseudoNuTildaEqn_
Definition: DASpalartAllmarasFv3.H:92
Foam::DASpalartAllmarasFv3::addModelResidualCon
virtual void addModelResidualCon(HashTable< List< List< word >>> &allCon) const
add the model residual connectivity to stateCon
Definition: DASpalartAllmarasFv3.C:335
Foam::DASpalartAllmarasFv3::calcLduResidualTurb
virtual void calcLduResidualTurb(volScalarField &nuTildaRes)
calculate the turbulence residual using LDU matrix
Definition: DASpalartAllmarasFv3.C:665
Foam::DASpalartAllmarasFv3::kappa_
dimensionedScalar kappa_
Definition: DASpalartAllmarasFv3.H:56
Foam::DASpalartAllmarasFv3::updateIntermediateVariables
virtual void updateIntermediateVariables()
update any intermediate variables that are dependent on state variables and are used in calcResiduals
Definition: DASpalartAllmarasFv3.C:274
Foam::DASpalartAllmarasFv3::fw
tmp< volScalarField > fw(const volScalarField &Stilda) const
Definition: DASpalartAllmarasFv3.C:184
Foam::DASpalartAllmarasFv3::constructPseudoNuTildaEqn
virtual void constructPseudoNuTildaEqn()
construct pseudoNuTildaEqn_, which is nuTildaEqn with swapped upper and lower arrays
Definition: DASpalartAllmarasFv3.C:595
Foam::DASpalartAllmarasFv3::sigmaNut_
dimensionedScalar sigmaNut_
Definition: DASpalartAllmarasFv3.H:55
Foam::DASpalartAllmarasFv3::invTranProdNuTildaEqn
virtual void invTranProdNuTildaEqn(const volScalarField &mySource, volScalarField &pseudoNuTilda)
Inverse transpose product, A_nuTilda^(-T)
Definition: DASpalartAllmarasFv3.C:520
Foam::DATurbulenceModel::U_
volVectorField & U_
velocity
Definition: DATurbulenceModel.H:86
Foam::DATurbulenceModel::nut_
volScalarField & nut_
turbulence viscosity
Definition: DATurbulenceModel.H:83
Foam::DASpalartAllmarasFv3::fv2
tmp< volScalarField > fv2(const volScalarField &chi, const volScalarField &fv1) const
Definition: DASpalartAllmarasFv3.C:164
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DASpalartAllmarasFv3::Cb2_
dimensionedScalar Cb2_
Definition: DASpalartAllmarasFv3.H:58
Foam::DATurbulenceModel::nu
tmp< volScalarField > nu() const
get the nu field
Definition: DATurbulenceModel.C:262
Foam::DATurbulenceModel::correctAlphat
void correctAlphat()
update alphat
Definition: DATurbulenceModel.C:206
solve
pseudoPEqn solve(solverDictP_)
Foam::DATurbulenceModel
Definition: DATurbulenceModel.H:61
Foam::DASpalartAllmarasFv3::chi
tmp< volScalarField > chi() const
Definition: DASpalartAllmarasFv3.C:152
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFvSource, 0)
Foam::DASpalartAllmarasFv3::pseudoNuTilda_
volScalarField pseudoNuTilda_
pseudoNuTilda_ and pseudoNuTildaEqn_ for solving adjoint equation
Definition: DASpalartAllmarasFv3.H:91
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::DASpalartAllmarasFv3::Cv2_
dimensionedScalar Cv2_
Definition: DASpalartAllmarasFv3.H:63
Foam::DASpalartAllmarasFv3::DASpalartAllmarasFv3
DASpalartAllmarasFv3(const word modelType, const fvMesh &mesh, const DAOption &daOption)
Definition: DASpalartAllmarasFv3.C:41
Foam::DASpalartAllmarasFv3::Cw2_
dimensionedScalar Cw2_
Definition: DASpalartAllmarasFv3.H:60
Foam::DASpalartAllmarasFv3::printInterval_
label printInterval_
time step interval to print residual
Definition: DASpalartAllmarasFv3.H:101
Foam::DASpalartAllmarasFv3::correctBoundaryConditions
virtual void correctBoundaryConditions()
update turbulence variable boundary values
Definition: DASpalartAllmarasFv3.C:263
Foam::DASpalartAllmarasFv3::DnuTildaEff
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
Definition: DASpalartAllmarasFv3.C:202
Foam::DASpalartAllmarasFv3::rhsSolvePseudoNuTildaEqn
virtual void rhsSolvePseudoNuTildaEqn(const volScalarField &nuTildaSource)
solve pseudoNuTildaEqn_ with overwritten rhs
Definition: DASpalartAllmarasFv3.C:633
Foam::DASpalartAllmarasFv3::correct
virtual void correct()
solve the residual equations and update the state
Definition: DASpalartAllmarasFv3.C:413
Foam::DATurbulenceModel::phase_
volScalarField phase_
phase field
Definition: DATurbulenceModel.H:92
Foam::DASpalartAllmarasFv3::correctStateResidualModelCon
virtual void correctStateResidualModelCon(List< List< word >> &stateCon) const
update the original variable connectivity for the adjoint state residuals in stateCon
Definition: DASpalartAllmarasFv3.C:285