adjointkOmegaSST.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2014-2023 PCOpt/NTUA
9  Copyright (C) 2014-2023 FOSS GP
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "adjointkOmegaSST.H"
30 #include "wallFvPatch.H"
33 #include "linear.H"
34 #include "reverseLinear.H"
37 #include "fvOptions.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 namespace incompressibleAdjoint
46 {
47 namespace adjointRASModels
48 {
49 
50 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
51 
54 
55 
56 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
57 
59 {
61  (
62  min
63  (
64  max
65  (
66  (scalar(1)/betaStar_)*sqrt(k())/(omega()*y_),
67  scalar(500)*nu()/(sqr(y_)*omega())
68  ),
70  ),
71  scalar(10)
72  );
73 
74  return tanh(pow4(arg1));
75 }
76 
77 
79 {
81  (
82  max
83  (
84  (scalar(2)/betaStar_)*sqrt(k())/(omega()*y_),
85  scalar(500)*nu()/(sqr(y_)*omega())
86  ),
87  scalar(100)
88  );
89 
90  return tanh(sqr(arg2));
91 }
92 
93 
95 (
96  const volScalarField& GbyNu0,
97  const volScalarField& F2,
98  const volScalarField& S2
99 ) const
100 {
101  return min
102  (
103  GbyNu0,
105  *max(a1_*omega(), b1_*F2*sqrt(S2))
106  );
107 }
108 
109 
111 (
112  const volScalarField::Internal& GbyNu0,
114  const volScalarField::Internal& S2
115 ) const
116 {
117  return min
118  (
119  GbyNu0,
120  (c1_/a1_)*betaStar_*omega()()
121  *max(a1_*omega()(), b1_*F2*sqrt(S2))
122  );
123 }
124 
125 
127 {
128  auto tzeroFirstCell =
130  (
131  IOobject
132  (
133  "zeroFirstCell",
134  runTime_.timeName(),
135  mesh_,
138  ),
139  mesh_,
141  );
142  volScalarField& zeroFirstCell = tzeroFirstCell.ref();
143 
145  label counter(0);
146  for (const fvPatchScalarField& omegab : omega().boundaryField())
147  {
148  const fvPatch& patch = omegab.patch();
149  if (isA<omegaWallFunctionFvPatchScalarField>(omegab))
150  {
151  const label patchi = patch.index();
152  const labelList& faceCells = patch.faceCells();
153  fvPatchScalarField& bf = zeroFirstCell.boundaryFieldRef()[patchi];
154  forAll(faceCells, faceI)
155  {
156  const label cellI = faceCells[faceI];
157 
158  zeroFirstCell[cellI] = 0.;
159  bf[faceI] = 0.;
160  firstCellIDs_[counter++] = cellI;
161  }
162  }
163  }
164  firstCellIDs_.setSize(counter);
166  zeroFirstCell.correctBoundaryConditions();
167 
168  return tzeroFirstCell;
169 }
170 
171 
173 {
174  const volVectorField& U = primalVars_.U();
175  const volVectorField& Ua = adjointVars_.UaInst();
176  word scheme("coeffsDiff");
177  tmp<volScalarField> tdR_dnut
178  (
179  dNutdbMult(U, Ua, nutRef(), scheme)
180  + dNutdbMult(k(), ka(), alphaK_, nutRef(), scheme)
183  );
184  volScalarField& dRdnut = tdR_dnut.ref();
185 
186  forAll(mesh_.boundary(), pI)
187  {
188  const fvPatch& patch = mesh_.boundary()[pI];
189  const fvPatchScalarField& nutb = nutRef().boundaryField()[pI];
190  const labelList& faceCells = patch.faceCells();
191  if (isA<nutkWallFunctionFvPatchScalarField>(nutb))
192  {
193  fvPatchScalarField& bf = dRdnut.boundaryFieldRef()[pI];
194  const scalarField kSnGrad(k().boundaryField()[pI].snGrad());
195  const scalarField omegaSnGrad
196  (
197  omega().boundaryField()[pI].snGrad()
198  );
199  const fvPatchScalarField& akb = alphaK_.boundaryField()[pI];
200  const fvPatchScalarField& aOmegab = alphaOmega_.boundaryField()[pI];
201  const vectorField USnGrad(U.boundaryField()[pI].snGrad());
202  const fvPatchTensorField& gradUb = gradU_.boundaryField()[pI];
203  const vectorField nf(mesh_.boundary()[pI].nf());
204  forAll(faceCells, fI)
205  {
206  const label cI(faceCells[fI]);
207  bf[fI] =
208  - wa()[cI]*zeroFirstCell_[cI]*aOmegab[fI]*omegaSnGrad[fI]
209  - ka()[cI]*akb[fI]*kSnGrad[fI]
210  - (Ua[cI] & (USnGrad[fI] + (dev2(gradUb[fI]) & nf[fI])));
211  }
212  }
213  }
214 
215  return tdR_dnut;
216 }
217 
218 
220 (
221  const volScalarField& F2,
222  const volScalarField& S,
223  const volScalarField& case_1_nut,
224  const volScalarField& case_2_nut,
225  const volScalarField& case_3_nut
226 ) const
227 {
228  return
229  (
230  - case_1_nut*k()/sqr(omega())
231  - a1_*k()/(b1_*S*F2*F2)*dF2_domega(F2, case_2_nut, case_3_nut)
232  );
233 }
234 
235 
237 (
238  const volScalarField& F2,
239  const volScalarField& S,
240  const volScalarField& case_2_nut
241 ) const
242 {
243  return
244  (
246  - a1_*k()/(b1_*S*F2*F2)*dF2_dk(F2, case_2_nut)
247  );
248 }
249 
250 
252 (
253  const volScalarField& F2,
254  const volScalarField& case_2_nut,
255  const volScalarField& case_3_nut
256 ) const
257 {
258  tmp<volScalarField> arg2 = min
259  (
260  max
261  (
262  (scalar(2)/betaStar_)*sqrt(k())/(omega()*y_),
263  scalar(500)*nu()/(sqr(y_)*omega())
264  ),
265  scalar(100)
266  );
267 
268  return
269  - scalar(2)*arg2*(1 - F2*F2)*
270  (
271  case_2_nut*scalar(2)*sqrt(k())/(betaStar_*sqr(omega())*y_)
272  + case_3_nut*scalar(500)*nu()/sqr(omega()*y_)
273  );
274 }
275 
276 
278 (
279  const volScalarField& F2,
280  const volScalarField& case_2_nut
281 ) const
282 {
283  tmp<volScalarField> arg2 = min
284  (
285  max
286  (
287  (scalar(2)/betaStar_)*sqrt(k())/(omega()*y_),
288  scalar(500)*nu()/(sqr(y_)*omega())
289  ),
290  scalar(100)
291  );
293  return
294  case_2_nut*scalar(2)*arg2*(1 - F2*F2)
295  /(betaStar_*omega()*y_*sqrt(k()));
296 }
297 
298 
300 {
301  return
302  (
304  *(
305  max(a1_*omega(), b1_*F2_*S_)
306  + case_1_nut_*a1_*omega()
307  + (scalar(1) - case_1_nut_)*omega()*b1_*S_
309  )
310  );
311 }
312 
313 
315 {
316  return
319 }
320 
321 
323 {
324  const volVectorField& U = primalVars_.U();
326  tmp<volScalarField> tGPrime = GbyNu(GbyNu0_, F2_, S2_);
328  word scheme("coeffsDiff");
329 
330  tmp<volScalarField> tdRdF1
331  (
332  nutRef()
333  *(
337  )
338  );
339  volScalarField& dRdF1 = tdRdF1.ref();
340 
341  typedef fixedValueFvPatchScalarField fixedValue;
342  typedef zeroGradientFvPatchScalarField zeroGrad;
343  typedef omegaWallFunctionFvPatchScalarField omegaWF;
344  forAll(mesh_.boundary(), pI)
345  {
346  const fvPatchScalarField& kb = k().boundaryField()[pI];
347  const fvPatchScalarField& omegab = omega().boundaryField()[pI];
348  fvPatchScalarField& dRdF1b = dRdF1.boundaryFieldRef()[pI];
349  if
350  (
351  isA<zeroGrad>(kb)
352  && (isA<zeroGrad>(omegab) || isA<omegaWF>(omegab))
353  )
354  {
355  dRdF1b = dRdF1b.patchInternalField();
356  }
357  else if (isA<fixedValue>(kb) && isA<fixedValue>(omegab)
358  )
359  {
360  // Note: might need revisiting
361  dRdF1b = dRdF1b.patchInternalField();
362  }
363  }
364 
365  volScalarField dR_dF1_coeffs
366  (
368  *(
369  (gamma1_ - gamma2_)*((2.0/3.0)*omega()*tdivU - tGPrime)
370  + omega()*omega()*(beta1_ - beta2_) + CDkOmega_
371  )
372  );
373 
374  forAll(mesh_.boundary(), pI)
375  {
376  const fvPatchScalarField& kb = k().boundaryField()[pI];
377  const fvPatchScalarField& omegab = omega().boundaryField()[pI];
378  fvPatchScalarField& dRdF1b = dR_dF1_coeffs.boundaryFieldRef()[pI];
379  if
380  (
381  isA<zeroGrad>(kb)
382  && (isA<zeroGrad>(omegab) || isA<omegaWF>(omegab))
383  )
384  {
385  dRdF1b = Zero;
386  }
387  else if (isA<fixedValue>(kb) && isA<fixedValue>(omegab))
388  {
389  dRdF1b = dRdF1b.patchInternalField();
390  }
391  }
393  dRdF1 += dR_dF1_coeffs;
394  return tdRdF1;
395 }
396 
397 
399 (
400  const volScalarField& arg1
401 ) const
402 {
403  return
404  (
405  4*pow3(arg1)*(scalar(1) - F1_*F1_)
406  *(
407  - case_1_F1_*sqrt(k())/(betaStar_*omega()*omega()*y_)
408  - case_2_F1_*scalar(500)*nu()/sqr(omega()*y_)
411  )
412  );
413 }
414 
415 
417 (
418  const volScalarField& arg1
419 ) const
420 {
421  return
422  (
423  - case_3_F1_*scalar(4)*pow3(arg1)*(scalar(1) - F1_*F1_)
424  *scalar(8)*k()*sqr(alphaOmega2_/(CDkOmegaPlus_*y_))/omega()*gradK_
425  );
426 }
427 
428 
430 {
431  tmp<volScalarField> arg1 = min
432  (
433  min
434  (
435  max
436  (
437  (scalar(1)/betaStar_)*sqrt(k())/(omega()*y_),
438  scalar(500)*nu()/(sqr(y_)*omega())
439  ),
441  ),
442  scalar(10)
443  );
444 
445  volScalarField dR_dF1(this->dR_dF1());
446  volScalarField dF1_domega(this->dF1_domega(arg1));
448 
449  surfaceScalarField dR_dGradOmega
450  (
451  interpolationScheme<vector>("div(dR_dGradOmega)")().
453  );
454 
455  return
456  (
458  - fvc::div(dR_dGradOmega)
459  );
460 }
461 
462 
464 {
465  tmp<volVectorField> tsource
466  (
468  );
469  volVectorField& source = tsource.ref();
470 
471  forAll(mesh_.boundary(), pI)
472  {
473  const fvPatchScalarField& omegab = omega().boundaryField()[pI];
474  fvPatchVectorField& sourceb = source.boundaryFieldRef()[pI];
475  if
476  (
477  isA<zeroGradientFvPatchScalarField>(omegab)
478  || isA<omegaWallFunctionFvPatchScalarField>(omegab)
479  )
480  {
481  sourceb = Zero;
482  }
483  else if (isA<fixedValueFvPatchScalarField>(omegab))
484  {
485  sourceb = sourceb.patchInternalField();
486  }
487  }
488 
489  surfaceScalarField sourceFaces
490  (
491  interpolationScheme<vector>("sourceAdjontkOmegaSST")().
492  interpolate(source) & mesh_.Sf()
493  );
494 
495  return
496  (
497  // Differentiation of omega in CDkOmega
499  // Differentiation of grad(omega) in CDkOmega
500  + fvc::div(sourceFaces)
501  );
502 }
503 
504 
506 {
507  tmp<volVectorField> tsource
508  (
510  );
511  volVectorField& source = tsource.ref();
512 
513  forAll(mesh_.boundary(), pI)
514  {
515  const fvPatchScalarField& kb = k().boundaryField()[pI];
516  fvPatchVectorField& sourceb = source.boundaryFieldRef()[pI];
517  if (isA<zeroGradientFvPatchScalarField>(kb))
518  {
519  sourceb = Zero;
520  }
521  else if (isA<fixedValueFvPatchScalarField>(kb))
522  {
523  sourceb = sourceb.patchInternalField();
524  }
525  }
526 
527  return
528  fvc::div
529  (
530  interpolationScheme<vector>("sourceAdjontkOmegaSST")().
531  interpolate(source) & mesh_.Sf()
532  );
533 }
534 
535 
537 (
538  const volScalarField& arg1
539 ) const
540 {
541  return
542  (
543  4*pow3(arg1)*(scalar(1) - F1_*F1_)
544  *(
545  case_1_F1_*0.5/(betaStar_*omega()*sqrt(k())*y_)
547  )
548  );
549 }
550 
551 
553 (
554  const volScalarField& arg1
555 ) const
556 {
557  return
558  (
559  - case_3_F1_*scalar(4)*pow3(arg1)*(scalar(1) - F1_*F1_)
560  *scalar(8)*k()*sqr(alphaOmega2_/(CDkOmegaPlus_*y_))/omega()
561  *gradOmega_
562  );
563 }
564 
565 
567 {
568  tmp<volScalarField> arg1 = min
569  (
570  min
571  (
572  max
573  (
574  (scalar(1)/betaStar_)*sqrt(k())/(omega()*y_),
575  scalar(500)*nu()/(sqr(y_)*omega())
576  ),
578  ),
579  scalar(10)
580  );
581 
582  volScalarField dR_dF1(this->dR_dF1());
583  volScalarField dF1_dk(this->dF1_dk(arg1));
584  volVectorField dF1_dGradK(this->dF1_dGradK(arg1));
585 
586  surfaceScalarField dR_dGradK
587  (
588  interpolationScheme<vector>("div(dR_dGradK)")().
590  );
591 
592  return
593  (
595  - fvc::div(dR_dGradK)
596  );
597 }
598 
599 
601 (
602  const volScalarField& primalField,
603  const volScalarField& adjointField,
604  const word& schemeName
605 ) const
606 {
608  (interpolationScheme<scalar>(schemeName));
609  surfaceScalarField snGradPrimal(fvc::snGrad(primalField)*mesh_.magSf());
610  surfaceScalarField flux(scheme().interpolate(adjointField)*snGradPrimal);
611 
612  forAll(mesh_.boundary(), pI)
613  {
614  const fvPatchScalarField& primalbf = primalField.boundaryField()[pI];
615  if (isA<zeroGradientFvPatchScalarField>(primalbf))
616  {
617  // Needless, but just to be safe
618  snGradPrimal.boundaryFieldRef()[pI] = Zero;
619  flux.boundaryFieldRef()[pI] = Zero;
620  }
621  else if (isA<fixedValueFvPatchScalarField>(primalbf))
622  {
623  // Note: to be potentially revisited
624  snGradPrimal.boundaryFieldRef()[pI] = Zero;
625  flux.boundaryFieldRef()[pI] = Zero;
626  }
627  }
628 
629  return (fvc::div(flux) - adjointField*fvc::div(snGradPrimal));
630 }
631 
632 
634 (
635  const volScalarField& primalField,
636  const volScalarField& adjointField,
637  const volScalarField& coeffField,
638  const volScalarField& bcField,
639  const word& schemeName
640 ) const
641 {
643  (interpolationScheme<scalar>(schemeName));
644  surfaceScalarField snGradPrimal(fvc::snGrad(primalField)*mesh_.magSf());
645  surfaceScalarField flux(scheme().interpolate(adjointField)*snGradPrimal);
646 
647  forAll(mesh_.boundary(), pI)
648  {
649  const fvPatchScalarField& bc = bcField.boundaryField()[pI];
650  if (isA<zeroGradientFvPatchScalarField>(bc))
651  {
652  const fvPatchScalarField& coeffFieldb =
653  coeffField.boundaryField()[pI];
654  snGradPrimal.boundaryFieldRef()[pI] *=
655  coeffFieldb/coeffFieldb.patchInternalField();
656  flux.boundaryFieldRef()[pI] = Zero;
657  }
658  else if (isA<fixedValueFvPatchScalarField>(bc))
659  {
660  snGradPrimal.boundaryFieldRef()[pI] = Zero;
661  flux.boundaryFieldRef()[pI] = Zero;
662  }
663  }
664 
665  return coeffField*(fvc::div(flux) - adjointField*fvc::div(snGradPrimal));
666 }
667 
668 
670 (
671  const volVectorField& U,
672  const volVectorField& Ua,
673  const volScalarField& nut,
674  const word& schemeName
675 ) const
676 {
678  (interpolationScheme<vector>(schemeName));
680  surfaceScalarField flux(scheme().interpolate(Ua) & snGradU);
681 
682  // Terms form the Laplacian part of the momentum stresses
683  forAll(mesh_.boundary(), pI)
684  {
685  const fvPatchScalarField& bc = nut.boundaryField()[pI];
686  if (isA<zeroGradientFvPatchScalarField>(bc))
687  {
688  flux.boundaryFieldRef()[pI] = Zero;
689  }
690  else if (isA<fixedValueFvPatchScalarField>(bc))
691  {
692  snGradU.boundaryFieldRef()[pI] = Zero;
693  flux.boundaryFieldRef()[pI] = Zero;
694  }
695  }
696 
697  // Terms form the tranpose gradient in the momentum stresses
698  const surfaceVectorField& Sf = mesh_.Sf();
699  surfaceTensorField fluxTranspose
700  (
701  reverseLinear<vector>(mesh_).interpolate(Ua)*Sf
702  );
703  forAll(mesh_.boundary(), pI)
704  {
705  const fvPatchVectorField& Ub = U.boundaryField()[pI];
706  if (!isA<coupledFvPatchVectorField>(Ub))
707  {
708  const vectorField Uai(Ua.boundaryField()[pI].patchInternalField());
709  const vectorField& Sfb = Sf.boundaryField()[pI];
710  fluxTranspose.boundaryFieldRef()[pI] = Uai*Sfb;
711  }
712  }
713  volScalarField M(dev2(gradU_) && fvc::div(fluxTranspose));
714  const DimensionedField<scalar, volMesh>& V = mesh_.V();
715  forAll(mesh_.boundary(), pI)
716  {
717  const fvPatchScalarField& bc = nut.boundaryField()[pI];
718  if (isA<zeroGradientFvPatchScalarField>(bc))
719  {
720  const vectorField Uai(Ua.boundaryField()[pI].patchInternalField());
721  const tensorField dev2GradU
722  (
723  dev2(gradU_.boundaryField()[pI].patchInternalField())
724  );
725  const vectorField& Sfb = Sf.boundaryField()[pI];
726  const labelList& faceCells = mesh_.boundary()[pI].faceCells();
727  forAll(faceCells, fI)
728  {
729  const label celli = faceCells[fI];
730  M[celli] -= ((Uai[fI] & dev2GradU[fI]) & Sfb[fI])/V[celli];
731  }
732  }
733  }
734  M.correctBoundaryConditions();
735 
736  //surfaceScalarField fluxTranspose =
737  // fvc::interpolate(UaGradU, schemeName) & mesh_.Sf();
738  //forAll(mesh_.boundary(), pI)
739  //{
740  // const fvPatchScalarField& bc = nut.boundaryField()[pI];
741  // const vectorField& Sf = mesh_.boundary()[pI].Sf();
742  // if (isA<zeroGradientFvPatchScalarField>(bc))
743  // {
744  // fluxTranspose.boundaryFieldRef()[pI] =
745  // (
746  // UaGradU.boundaryField()[pI].patchInternalField()
747  // - (
748  // Ua.boundaryField()[pI].patchInternalField()
749  // & gradU_.boundaryField()[pI]
750  // )
751  // ) & Sf;
752  // }
753  // else if (isA<fixedValueFvPatchScalarField>(bc))
754  // {
755  // fluxTranspose.boundaryFieldRef()[pI] =
756  // UaGradU.boundaryField()[pI].patchInternalField() & Sf;
757  // }
758  //}
759  return
760  fvc::div(flux) - (Ua & fvc::div(snGradU)) + M;
761  //fvc::div(flux + fluxTranspose) - (Ua & fvc::div(snGradU));
762 }
763 
764 
766 (
767  const volScalarField& primalField,
768  const volScalarField& adjointField
769 ) const
770 {
771  // Grab the interpolation scheme from the primal convection term
772  tmp<surfaceInterpolationScheme<scalar>> primalInterpolationScheme
773  (
774  convectionScheme(primalField.name())
775  );
776 
777  surfaceVectorField interpolatedPrimal
778  (
779  primalInterpolationScheme().interpolate(primalField)*mesh_.Sf()
780  );
782  (
783  //reverseLinear<scalar>(mesh_).interpolate(adjointField)
784  linear<scalar>(mesh_).interpolate(adjointField)
785  *interpolatedPrimal
786  );
787 
788  const volVectorField& U = primalVars_.U();
789  forAll(mesh_.boundary(), pI)
790  {
791  const fvPatchVectorField& bc = U.boundaryField()[pI];
792  if (isA<zeroGradientFvPatchVectorField>(bc))
793  {
794  flux.boundaryFieldRef()[pI] = Zero;
795  }
796  else if (isA<fixedValueFvPatchVectorField>(bc))
797  {
798  interpolatedPrimal.boundaryFieldRef()[pI] = Zero;
799  flux.boundaryFieldRef()[pI] = Zero;
800  }
801  }
802 
803  return (-fvc::div(flux) + adjointField*fvc::div(interpolatedPrimal));
804 }
805 
806 
808 (
809  tmp<volSymmTensorField>& GbyNuMult
810 ) const
811 {
813  (
815  );
816 
817  const volVectorField& U = primalVars_.U();
818  forAll(mesh_.boundary(), pI)
819  {
820  const fvPatchVectorField& bc = U.boundaryField()[pI];
821  if (isA<zeroGradientFvPatchVectorField>(bc))
822  {
823  flux.boundaryFieldRef()[pI] = Zero;
824  }
825  else if (isA<fixedValueFvPatchVectorField>(bc))
826  {
827  flux.boundaryFieldRef()[pI] =
828  mesh_.boundary()[pI].Sf()
829  & GbyNuMult().boundaryField()[pI].patchInternalField();
830  }
831  }
832 
833  return fvc::div(flux);
834 }
835 
836 
838 (
839  tmp<volScalarField>& divUMult
840 ) const
841 {
843  (
845  );
846 
847  const volVectorField& U = primalVars_.U();
848  forAll(mesh_.boundary(), pI)
849  {
850  const fvPatchVectorField& bc = U.boundaryField()[pI];
851  if (isA<zeroGradientFvPatchVectorField>(bc))
852  {
853  flux.boundaryFieldRef()[pI] = Zero;
854  }
855  else if (isA<fixedValueFvPatchVectorField>(bc))
856  {
857  flux.boundaryFieldRef()[pI] =
858  mesh_.boundary()[pI].Sf()
859  *divUMult().boundaryField()[pI].patchInternalField();
860  }
861  }
862 
863  divUMult.clear();
864 
865  return -fvc::div(flux);
866 }
867 
868 
870 (
871  const volScalarField& primalField,
872  const volScalarField& adjointField,
873  const volScalarField& coeffField
874 ) const
875 {
876  // Note: we could grab the snGrad scheme from the diffusion term directly
877  surfaceScalarField snGradPrimal(fvc::snGrad(primalField)*mesh_.magSf());
879  (
880  reverseLinear<scalar>(mesh_).interpolate(adjointField)*snGradPrimal
881  );
882 
883  const volVectorField& U = primalVars_.U();
884  forAll(mesh_.boundary(), pI)
885  {
886  const fvPatchVectorField& bc = U.boundaryField()[pI];
887  if (!isA<coupledFvPatchVectorField>(bc))
888  {
889  flux.boundaryFieldRef()[pI] = Zero;
890  snGradPrimal.boundaryFieldRef()[pI] = Zero;
891  }
892  }
893  return (fvc::div(flux) - adjointField*fvc::div(snGradPrimal))*coeffField;
894 }
895 
896 
898 (
899  tmp<volScalarField>& mult,
900  const volScalarField& F2,
901  const volScalarField& S,
902  const volScalarField& case_1_nut,
903  const volTensorField& gradU
904 ) const
905 {
907  (
908  mult*a1_*k()*(1 - case_1_nut)/(b1_*F2*S*S*S)*twoSymm(gradU)
909  );
910  M.correctBoundaryConditions();
911  mult.clear();
912 
913  surfaceVectorField returnFieldFlux
914  (
916  );
917 
918  const volVectorField& U = primalVars_.U();
919  forAll(mesh_.boundary(), pI)
920  {
921  const fvPatchVectorField& bc = U.boundaryField()[pI];
922  if (isA<zeroGradientFvPatchVectorField>(bc))
923  {
924  returnFieldFlux.boundaryFieldRef()[pI] = Zero;
925  }
926  else if (isA<fixedValueFvPatchVectorField>(bc))
927  {
928  returnFieldFlux.boundaryFieldRef()[pI] =
929  mesh_.boundary()[pI].Sf()
930  & M.boundaryField()[pI].patchInternalField();
931  }
932  }
933 
934  // Note: potentially missing contributions form patches with a calculated
935  // nut bc
936 
937  return fvc::div(returnFieldFlux);
938 }
939 
940 
942 (
943  fvScalarMatrix& kaEqn,
944  const volScalarField& dR_dnut
945 )
946 {
947  // Add source term from the differentiation of nutWallFunction
948  scalarField& source = kaEqn.source();
950  const volScalarField& ka = this->ka();
951  const volScalarField& wa = this->wa();
952  const volScalarField& k = this->k();
953  const volScalarField& omega = this->omega();
954  forAll(nutRef().boundaryFieldRef(), patchi)
955  {
956  fvPatchScalarField& nutWall = nutRef().boundaryFieldRef()[patchi];
957  if (isA<nutkWallFunctionFvPatchScalarField>(nutWall))
958  {
959  const fvPatch& patch = mesh_.boundary()[patchi];
960  const scalarField& magSf = patch.magSf();
961 
964 
965  const scalarField& y = turbModel().y()[patchi];
966  const tmp<scalarField> tnuw = turbModel().nu(patchi);
967  const scalarField& nuw = tnuw();
968 
969  const nutWallFunctionFvPatchScalarField& nutWF =
970  refCast<nutWallFunctionFvPatchScalarField>(nutWall);
971  const wallFunctionCoefficients& wallCoeffs = nutWF.wallCoeffs();
972  const scalar Cmu = wallCoeffs.Cmu();
973  const scalar kappa = wallCoeffs.kappa();
974  const scalar E = wallCoeffs.E();
975  const scalar yPlusLam = wallCoeffs.yPlusLam();
976 
977  const scalar Cmu25 = pow025(Cmu);
978 
979  const labelList& faceCells = patch.faceCells();
980  const fvPatchScalarField& dR_dnutw =
981  dR_dnut.boundaryField()[patchi];
982  const fvPatchScalarField& omegaw = omega.boundaryField()[patchi];
983  bool addTermsFromOmegaWallFuction
984  (isA<omegaWallFunctionFvPatchScalarField>(omegaw));
985 
986  const fvPatchVectorField& Uw =
987  primalVars_.U().boundaryField()[patchi];
988  const scalarField magGradUw(mag(Uw.snGrad()));
989  forAll(nuw, facei)
990  {
991  const label celli = faceCells[facei];
992 
993  const scalar sqrtkCell(sqrt(k[celli]));
994  const scalar yPlus = Cmu25*y[facei]*sqrtkCell/nuw[facei];
995  const scalar logEyPlus = log(E*yPlus);
996  const scalar dnut_dyPlus =
997  nuw[facei]*kappa*(logEyPlus - 1)/sqr(logEyPlus);
998  const scalar dyPlus_dk =
999  Cmu25*y[facei]/(2*nuw[facei]*sqrtkCell);
1000  const scalar dnut_dk = dnut_dyPlus*dyPlus_dk;
1001 
1002  if (yPlusLam < yPlus)
1003  {
1004  // Terms from nutkWallFunction
1005  source[celli] -= dR_dnutw[facei]*dnut_dk*magSf[facei];
1006  }
1007  if (addTermsFromOmegaWallFuction)
1008  {
1009  const scalar denom(Cmu25*kappa*y[facei]);
1010  const scalar omegaLog(sqrtkCell/denom);
1011  // Terms from omegaLog in omegaWallFunction
1012  source[celli] +=
1013  wa[celli]*omegaLog/omega[celli]
1014  /(2*sqrtkCell*denom);
1015 
1016  // Terms from G at the first cell off the wall
1017  source[celli] +=
1018  case_1_Pk_[celli]*ka[celli]*V[celli]
1019  *(
1020  (nutWall[facei] + nuw[facei])
1021  *magGradUw[facei]
1022  *Cmu25
1023  /(2.0*sqrtkCell*kappa*y[facei])
1024  );
1025 
1026  if (yPlusLam < yPlus)
1027  {
1028  source[celli] +=
1029  case_1_Pk_[celli]*ka[celli]*V[celli]
1030  *dnut_dk
1031  *magGradUw[facei]
1032  *Cmu25*sqrtkCell
1033  /(kappa*y[facei]);
1034  }
1035  }
1036  }
1037  }
1038  }
1039 }
1040 
1041 
1043 {
1045  {
1046  Info<< "Updating primal-based fields of the adjoint turbulence "
1047  << "model ..." << endl;
1048 
1049  const volVectorField& U = primalVars_.U();
1050  gradU_ = fvc::grad(U);
1051  gradOmega_ = fvc::grad(omega());
1052  gradK_ = fvc::grad(k());
1053 
1054  S2_ = 2*magSqr(symm(gradU_))
1056  S_ = sqrt(S2_);
1058 
1059  // Instead of computing G directly here, delegate to RASModelVariables
1060  // to make sure G is computed consistently when the primal fields are
1061  // averaged. The local value computed here is just used to update
1062  // the switch fields
1063  /*
1064  volScalarField G(primalVars_.turbulence()->GName(), nutRef()*GbyNu0_);
1065  omega().correctBoundaryConditions();
1066  */
1068  (
1069  IOobject
1070  (
1071  IOobject::scopedName(type(), "G"),
1072  mesh_.time().timeName(),
1073  mesh_,
1076  ),
1077  mesh_,
1079  );
1080  G.ref() = primalVars_.RASModelVariables()->G()();
1081 
1082  CDkOmega_ =
1083  (2*alphaOmega2_)*(gradK_ & gradOmega_)/omega();
1084  CDkOmegaPlus_ = max
1085  (
1086  CDkOmega_,
1087  dimensionedScalar("1.0e-10", dimless/sqr(dimTime), 1.0e-10)
1088  );
1089  F1_ = F1();
1090  F2_ = F2();
1091  case_1_Pk_ = neg(G - c1_*betaStar_*k()*omega());
1092  case_2_Pk_ = pos0(G - c1_*betaStar_*k()*omega());
1093  case_3_Pk_ = neg(a1_*omega() - b1_*F2_*S_);
1094 
1095 
1096  alphaK_ = alphaK(F1_);
1098  beta_ = beta(F1_);
1099  gamma_ = gamma(F1_);
1100 
1101  // Switch fields for F1
1102  {
1103  volScalarField arg1_C3
1104  (
1105  (scalar(1)/betaStar_)*sqrt(k())/(omega()*y_)
1106  - scalar(500)*nu()/(sqr(y_)*omega())
1107  );
1108  volScalarField arg1_C2
1109  (
1110  max
1111  (
1112  (scalar(1)/betaStar_)*sqrt(k())/(omega()*y_),
1113  scalar(500)*nu()/(sqr(y_)*omega())
1114  )
1115  - (4*alphaOmega2_)*k()/(CDkOmegaPlus_*sqr(y_))
1116  );
1117  volScalarField arg1_C1
1118  (
1119  min
1120  (
1121  max
1122  (
1123  (scalar(1)/betaStar_)*sqrt(k())/(omega()*y_),
1124  scalar(500)*nu()/(sqr(y_)*omega())
1125  ),
1126  (4*alphaOmega2_)*k()/(CDkOmegaPlus_*sqr(y_))
1127  ) - scalar(10)
1128  );
1129  volScalarField CDkOmegaPlus_C1
1130  (
1131  CDkOmega_
1132  - dimensionedScalar("1.0e-10", dimless/sqr(dimTime), 1.0e-10)
1133  );
1134 
1135  case_1_F1_ = pos(arg1_C3)*neg(arg1_C2)*neg(arg1_C1);
1136  case_2_F1_ = neg0(arg1_C3)*neg(arg1_C2)*neg(arg1_C1);
1137  case_3_F1_ = pos0(arg1_C2)*neg(arg1_C1)*pos(CDkOmegaPlus_C1);
1138  case_4_F1_ = pos0(arg1_C2)*neg(arg1_C1);
1139  }
1140 
1141  // Switch fields for nut
1142  {
1143  volScalarField nut_C1(a1_*omega() - b1_*F2_*S_);
1144  volScalarField arg2_C2
1145  (
1146  (scalar(2)/betaStar_)*sqrt(k())/(omega()*y_)
1147  - scalar(500)*nu()/(sqr(y_)*omega())
1148  );
1149  volScalarField arg2_C1
1150  (
1151  max
1152  (
1153  (scalar(2)/betaStar_)*sqrt(k())/(omega()*y_),
1154  scalar(500)*nu()/(sqr(y_)*omega())
1155  ) - scalar(100)
1156  );
1157 
1158  case_1_nut_ = pos(nut_C1);
1159  case_2_nut_ = neg0(nut_C1)*pos(arg2_C2)*neg(arg2_C1);
1160  case_3_nut_ = neg0(nut_C1)*neg0(arg2_C2)*neg(arg2_C1);
1161  }
1162 
1163  {
1164  volScalarField GPrime_C1
1165  (
1166  GbyNu0_
1167  - (c1_/a1_)*betaStar_*omega()*max(a1_*omega(), b1_*F2_*S_)
1168  );
1169  case_1_GPrime_ = neg(GPrime_C1);
1170  case_2_GPrime_ = pos0(GPrime_C1);
1171  }
1172 
1173  dnut_domega_ =
1177  DkEff_ = DkEff(F1_);
1179  changedPrimalSolution_ = false;
1180  }
1181 }
1182 
1183 
1185 (
1186  const word& varName
1187 ) const
1188 {
1190  const surfaceScalarField& phiInst = primalVars_.phiInst();
1191  word divEntry("div(" + phiInst.name() + ',' + varName +')');
1192  ITstream& divScheme = mesh_.divScheme(divEntry);
1193  // Skip the first entry which might be 'bounded' or 'Gauss'.
1194  // If it is 'bounded', skip the second entry as well
1195  word discarded(divScheme);
1196  if (discarded == "bounded")
1197  {
1198  discarded = word(divScheme);
1199  }
1201 }
1202 
1203 
1204 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1205 
1206 adjointkOmegaSST::adjointkOmegaSST
1207 (
1208  incompressibleVars& primalVars,
1209  incompressibleAdjointMeanFlowVars& adjointVars,
1210  objectiveManager& objManager,
1211  const word& adjointTurbulenceModelName,
1212  const word& modelName
1213 )
1214 :
1216  (
1217  modelName,
1218  primalVars,
1219  adjointVars,
1220  objManager,
1221  adjointTurbulenceModelName
1222  ),
1223 
1224  kappa_
1225  (
1226  dimensioned<scalar>::getOrAddToDict
1227  (
1228  "kappa",
1229  coeffDict_,
1230  0.41
1231  )
1232  ),
1233  alphaK1_
1234  (
1235  dimensioned<scalar>::getOrAddToDict
1236  (
1237  "alphaK1",
1238  this->coeffDict_,
1239  0.85
1240  )
1241  ),
1242  alphaK2_
1243  (
1244  dimensioned<scalar>::getOrAddToDict
1245  (
1246  "alphaK2",
1247  this->coeffDict_,
1248  1.0
1249  )
1250  ),
1251  alphaOmega1_
1252  (
1253  dimensioned<scalar>::getOrAddToDict
1254  (
1255  "alphaOmega1",
1256  this->coeffDict_,
1257  0.5
1258  )
1259  ),
1260  alphaOmega2_
1261  (
1262  dimensioned<scalar>::getOrAddToDict
1263  (
1264  "alphaOmega2",
1265  this->coeffDict_,
1266  0.856
1267  )
1268  ),
1269  gamma1_
1270  (
1271  dimensioned<scalar>::getOrAddToDict
1272  (
1273  "gamma1",
1274  this->coeffDict_,
1275  5.0/9.0
1276  )
1277  ),
1278  gamma2_
1279  (
1280  dimensioned<scalar>::getOrAddToDict
1281  (
1282  "gamma2",
1283  this->coeffDict_,
1284  0.44
1285  )
1286  ),
1287  beta1_
1288  (
1289  dimensioned<scalar>::getOrAddToDict
1290  (
1291  "beta1",
1292  this->coeffDict_,
1293  0.075
1294  )
1295  ),
1296  beta2_
1297  (
1298  dimensioned<scalar>::getOrAddToDict
1299  (
1300  "beta2",
1301  this->coeffDict_,
1302  0.0828
1303  )
1304  ),
1305  betaStar_
1306  (
1307  dimensioned<scalar>::getOrAddToDict
1308  (
1309  "betaStar",
1310  this->coeffDict_,
1311  0.09
1312  )
1313  ),
1314  a1_
1315  (
1316  dimensioned<scalar>::getOrAddToDict
1317  (
1318  "a1",
1319  this->coeffDict_,
1320  0.31
1321  )
1322  ),
1323  b1_
1324  (
1325  dimensioned<scalar>::getOrAddToDict
1326  (
1327  "b1",
1328  this->coeffDict_,
1329  1.0
1330  )
1331  ),
1332  c1_
1333  (
1334  dimensioned<scalar>::getOrAddToDict
1335  (
1336  "c1",
1337  this->coeffDict_,
1338  10.0
1339  )
1340  ),
1341  F3_
1342  (
1343  Switch::getOrAddToDict
1344  (
1345  "F3",
1346  this->coeffDict_,
1347  false
1348  )
1349  ),
1350 
1351  y_(primalVars_.RASModelVariables()().d()),
1352 
1353  //Primal Gradient Fields
1354  gradU_
1355  (
1356  IOobject
1357  (
1358  "rasModel::gradU",
1359  runTime_.timeName(),
1360  mesh_,
1361  IOobject::NO_READ,
1362  IOobject::NO_WRITE
1363  ),
1364  mesh_,
1366  ),
1367  gradOmega_
1368  (
1369  IOobject
1370  (
1371  "rasModel::gradOmega",
1372  runTime_.timeName(),
1373  mesh_,
1374  IOobject::NO_READ,
1375  IOobject::NO_WRITE
1376  ),
1377  mesh_,
1379  ),
1380  gradK_
1381  (
1382  IOobject
1383  (
1384  "rasModel::gradK",
1385  runTime_.timeName(),
1386  mesh_,
1387  IOobject::NO_READ,
1388  IOobject::NO_WRITE
1389  ),
1390  mesh_,
1392  ),
1393 
1394  S2_
1395  (
1396  IOobject
1397  (
1398  "S2",
1399  runTime_.timeName(),
1400  mesh_,
1401  IOobject::NO_READ,
1402  IOobject::NO_WRITE
1403  ),
1404  mesh_,
1406  ),
1407  S_
1408  (
1409  IOobject
1410  (
1411  "kOmegaSST_S",
1412  runTime_.timeName(),
1413  mesh_,
1414  IOobject::NO_READ,
1415  IOobject::NO_WRITE
1416  ),
1417  mesh_,
1419  ),
1420  GbyNu0_
1421  (
1422  IOobject
1423  (
1424  "adjointRASModel::GbyNu0",
1425  runTime_.timeName(),
1426  mesh_,
1427  IOobject::NO_READ,
1428  IOobject::NO_WRITE
1429  ),
1430  mesh_,
1432  ),
1433  CDkOmega_
1434  (
1435  IOobject
1436  (
1437  "CDkOmega_",
1438  runTime_.timeName(),
1439  mesh_,
1440  IOobject::NO_READ,
1441  IOobject::NO_WRITE
1442  ),
1443  mesh_,
1445  ),
1446  CDkOmegaPlus_
1447  (
1448  IOobject
1449  (
1450  "CDkOmegaPlus",
1451  runTime_.timeName(),
1452  mesh_,
1453  IOobject::NO_READ,
1454  IOobject::NO_WRITE
1455  ),
1456  mesh_,
1458  ),
1459  F1_
1460  (
1461  IOobject
1462  (
1463  "F1",
1464  runTime_.timeName(),
1465  mesh_,
1466  IOobject::NO_READ,
1467  IOobject::NO_WRITE
1468  ),
1469  mesh_,
1471  ),
1472  F2_
1473  (
1474  IOobject
1475  (
1476  "F2",
1477  runTime_.timeName(),
1478  mesh_,
1479  IOobject::NO_READ,
1480  IOobject::NO_WRITE
1481  ),
1482  mesh_,
1484  ),
1485  // Model Field coefficients
1486  alphaK_
1487  (
1488  IOobject
1489  (
1490  "alphaK",
1491  runTime_.timeName(),
1492  mesh_,
1493  IOobject::NO_READ,
1494  IOobject::NO_WRITE
1495  ),
1496  mesh_,
1498  ),
1499  alphaOmega_
1500  (
1501  IOobject
1502  (
1503  "alphaOmega",
1504  runTime_.timeName(),
1505  mesh_,
1506  IOobject::NO_READ,
1507  IOobject::NO_WRITE
1508  ),
1509  mesh_,
1511  ),
1512  beta_
1513  (
1514  IOobject
1515  (
1516  "beta",
1517  runTime_.timeName(),
1518  mesh_,
1519  IOobject::NO_READ,
1520  IOobject::NO_WRITE
1521  ),
1522  mesh_,
1524  ),
1525  gamma_
1526  (
1527  IOobject
1528  (
1529  "gamma",
1530  runTime_.timeName(),
1531  mesh_,
1532  IOobject::NO_READ,
1533  IOobject::NO_WRITE
1534  ),
1535  mesh_,
1537  ),
1538 
1539  case_1_F1_
1540  (
1541  IOobject
1542  (
1543  "case_1_F1",
1544  runTime_.timeName(),
1545  mesh_,
1546  IOobject::NO_READ,
1547  IOobject::NO_WRITE
1548  ),
1549  mesh_,
1551  ),
1552  case_2_F1_
1553  (
1554  IOobject
1555  (
1556  "case_2_F1",
1557  runTime_.timeName(),
1558  mesh_,
1559  IOobject::NO_READ,
1560  IOobject::NO_WRITE
1561  ),
1562  mesh_,
1564  ),
1565  case_3_F1_
1566  (
1567  IOobject
1568  (
1569  "case_3_F1",
1570  runTime_.timeName(),
1571  mesh_,
1572  IOobject::NO_READ,
1573  IOobject::NO_WRITE
1574  ),
1575  mesh_,
1577  ),
1578  case_4_F1_
1579  (
1580  IOobject
1581  (
1582  "case_4_F1",
1583  runTime_.timeName(),
1584  mesh_,
1585  IOobject::NO_READ,
1586  IOobject::NO_WRITE
1587  ),
1588  mesh_,
1590  ),
1591  case_1_Pk_
1592  (
1593  IOobject
1594  (
1595  "case_1_Pk",
1596  runTime_.timeName(),
1597  mesh_,
1598  IOobject::NO_READ,
1599  IOobject::NO_WRITE
1600  ),
1601  mesh_,
1603  ),
1604  case_2_Pk_
1605  (
1606  IOobject
1607  (
1608  "case_2_Pk",
1609  runTime_.timeName(),
1610  mesh_,
1611  IOobject::NO_READ,
1612  IOobject::NO_WRITE
1613  ),
1614  mesh_,
1616  ),
1617  case_3_Pk_
1618  (
1619  IOobject
1620  (
1621  "case_3_Pk",
1622  runTime_.timeName(),
1623  mesh_,
1624  IOobject::NO_READ,
1625  IOobject::NO_WRITE
1626  ),
1627  mesh_,
1629  ),
1630 
1631  case_1_nut_
1632  (
1633  IOobject
1634  (
1635  "case_1_nut",
1636  runTime_.timeName(),
1637  mesh_,
1638  IOobject::NO_READ,
1639  IOobject::NO_WRITE
1640  ),
1641  mesh_,
1643  ),
1644  case_2_nut_
1645  (
1646  IOobject
1647  (
1648  "case_2_nut",
1649  runTime_.timeName(),
1650  mesh_,
1651  IOobject::NO_READ,
1652  IOobject::NO_WRITE
1653  ),
1654  mesh_,
1656  ),
1657  case_3_nut_
1658  (
1659  IOobject
1660  (
1661  "case_3_nut",
1662  runTime_.timeName(),
1663  mesh_,
1664  IOobject::NO_READ,
1665  IOobject::NO_WRITE
1666  ),
1667  mesh_,
1669  ),
1670  case_1_GPrime_
1671  (
1672  IOobject
1673  (
1674  "case_1_GPrime",
1675  runTime_.timeName(),
1676  mesh_,
1677  IOobject::NO_READ,
1678  IOobject::NO_WRITE
1679  ),
1680  mesh_,
1682  ),
1683  case_2_GPrime_
1684  (
1685  IOobject
1686  (
1687  "case_2_GPrime",
1688  runTime_.timeName(),
1689  mesh_,
1690  IOobject::NO_READ,
1691  IOobject::NO_WRITE
1692  ),
1693  mesh_,
1695  ),
1696 
1697  // Zero 1rst cell field
1698  firstCellIDs_(0),
1699  zeroFirstCell_(zeroFirstCell()),
1700 
1701  // Turbulence model multipliers
1702  dnut_domega_
1703  (
1704  IOobject
1705  (
1706  "dnut_domega",
1707  runTime_.timeName(),
1708  mesh_,
1709  IOobject::NO_READ,
1710  IOobject::NO_WRITE
1711  ),
1712  mesh_,
1714  ),
1715  dnut_dk_
1716  (
1717  IOobject
1718  (
1719  "dnut_dk",
1720  runTime_.timeName(),
1721  mesh_,
1722  IOobject::NO_READ,
1723  IOobject::NO_WRITE
1724  ),
1725  mesh_,
1727  ),
1728  DOmegaEff_
1729  (
1730  IOobject
1731  (
1732  "DomegaEff",
1733  runTime_.timeName(),
1734  mesh_,
1735  IOobject::NO_READ,
1736  IOobject::NO_WRITE
1737  ),
1738  mesh_,
1739  dimensionedScalar(nutRef().dimensions(), Zero)
1740  ),
1741  DkEff_
1742  (
1743  IOobject
1744  (
1745  "DkEff",
1746  runTime_.timeName(),
1747  mesh_,
1748  IOobject::NO_READ,
1749  IOobject::NO_WRITE
1750  ),
1751  mesh_,
1752  dimensionedScalar(nutRef().dimensions(), Zero)
1753  )
1754 {
1756  adjointTMVariablesBaseNames_[0] = "ka";
1757  adjointTMVariablesBaseNames_[1] = "wa";
1758  // Read in adjoint fields
1760  (
1762  mesh_,
1763  "ka",
1764  adjointVars.solverName(),
1765  adjointVars.useSolverNameForFields()
1766  );
1768  (
1770  mesh_,
1771  "wa",
1772  adjointVars.solverName(),
1773  adjointVars.useSolverNameForFields()
1774  );
1775 
1776  setMeanFields();
1777 
1778  // No sensitivity contributions from the adjoint to the eikonal equation
1779  // for the moment
1780  includeDistance_ = false;
1781 
1782  // Update the primal related fields here so that functions computing
1783  // sensitivities have the updated fields in case of continuation
1785 }
1786 
1787 
1788 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1789 
1792  const volVectorField& Ua = adjointVars_.UaInst();
1793  return devReff(Ua);
1794 }
1795 
1796 
1798 (
1799  const volVectorField& U
1800 ) const
1801 {
1802  return
1804  (
1805  IOobject
1806  (
1807  "devRhoReff",
1808  runTime_.timeName(),
1809  mesh_,
1812  ),
1814  );
1815 }
1816 
1817 
1819 {
1820  tmp<volScalarField> tnuEff = nuEff();
1821  const volScalarField& nuEff = tnuEff();
1822 
1823  return
1824  (
1825  - fvm::laplacian(nuEff, Ua)
1826  - fvc::div(nuEff*dev(fvc::grad(Ua)().T()))
1827  );
1828 
1829  /* WIP
1830  const volVectorField& U = primalVars_.U();
1831  const surfaceVectorField& Sf = mesh_.Sf();
1832  tmp<surfaceTensorField> tflux =
1833  reverseLinear<vector>(mesh_).interpolate(Ua)*Sf;
1834  surfaceTensorField& flux = tflux.ref();
1835  forAll(mesh_.boundary(), pI)
1836  {
1837  const fvPatchVectorField& Ub = U.boundaryField()[pI];
1838  if (!isA<coupledFvPatchVectorField>(Ub))
1839  {
1840  const vectorField Uai = Ua.boundaryField()[pI].patchInternalField();
1841  const vectorField& Sfb = Sf.boundaryField()[pI];
1842  flux.boundaryFieldRef()[pI] = Uai*Sfb;
1843  }
1844  }
1845  volTensorField M(nuEff*dev2(fvc::div(flux)));
1846  const DimensionedField<scalar, volMesh>& V = mesh_.V();
1847 
1848  forAll(mesh_.boundary(), pI)
1849  {
1850  const fvPatchVectorField& Ub = U.boundaryField()[pI];
1851  if (!isA<coupledFvPatchVectorField>(Ub))
1852  {
1853  const fvPatchScalarField& nuEffb = nuEff.boundaryField()[pI];
1854  const vectorField nf = mesh_.boundary()[pI].nf();
1855  const vectorField Uai = Ua.boundaryField()[pI].patchInternalField();
1856  const labelList& faceCells = mesh_.boundary()[pI].faceCells();
1857  const vectorField& Sfb = Sf.boundaryField()[pI];
1858 
1859  forAll(faceCells, fI)
1860  {
1861  const label celli = faceCells[fI];
1862  const tensor t(dev2(Uai[fI]*Sfb[fI]));
1863  M[celli] -= nuEffb[fI]*(t - nf[fI]*(nf[fI] & t))/V[celli];
1864  }
1865  }
1866  }
1867  M.correctBoundaryConditions();
1868 
1869  surfaceVectorField returnFlux =
1870  - (Sf & reverseLinear<tensor>(mesh_).interpolate(M));
1871  forAll(mesh_.boundary(), pI)
1872  {
1873  const fvPatchVectorField& Ub = U.boundaryField()[pI];
1874  if (isA<zeroGradientFvPatchVectorField>(Ub))
1875  {
1876  returnFlux.boundaryFieldRef()[pI] = Zero;
1877  }
1878  else if (isA<fixedValueFvPatchVectorField>(Ub))
1879  {
1880  const scalarField& deltaCoeffs = mesh_.boundary()[pI].deltaCoeffs();
1881  const fvPatchScalarField& nuEffb = nuEff.boundaryField()[pI];
1882  const vectorField Uai = Ua.boundaryField()[pI].patchInternalField();
1883  const vectorField nf = mesh_.boundary()[pI].nf();
1884  const vectorField& Sfb = Sf.boundaryField()[pI];
1885 
1886  returnFlux.boundaryFieldRef()[pI] =
1887  - (Sfb & M.boundaryField()[pI].patchInternalField())
1888  + nuEffb*deltaCoeffs*(nf & dev2(Uai*Sfb));
1889  }
1890  }
1891 
1892  return
1893  (
1894  - fvm::laplacian(nuEff, Ua)
1895  + fvc::div(returnFlux)
1896  );
1897  */
1898 }
1899 
1902 {
1903  return (ka()*gradK_ + wa()*gradOmega_);
1904 }
1905 
1906 
1908 {
1909  tmp<volVectorField> tmeanFlowSource
1910  (
1912  (
1913  IOobject
1914  (
1915  "adjointMeanFlowSource" + type(),
1916  mesh_.time().timeName(),
1917  mesh_,
1920  ),
1921  mesh_,
1923  )
1924  );
1925  volVectorField& meanFlowSource = tmeanFlowSource.ref();
1926 
1927  // Contributions from the convection terms of the turbulence model
1928  meanFlowSource +=
1930  + convectionMeanFlowSource(k(), ka());
1931 
1932  // Contributions from GbyNu, including gradU
1933  tmp<volSymmTensorField> twoSymmGradU(twoSymm(gradU_));
1934  tmp<volSymmTensorField> GbyNuMult
1935  (
1936  // First part of GPrime and G from Pk
1937  2.*dev(twoSymmGradU())*zeroFirstCell_
1939  // Second part of GPrime
1940  + twoSymmGradU()*wa()*zeroFirstCell_*gamma_*case_2_GPrime_
1942  );
1943  twoSymmGradU.clear();
1944  meanFlowSource += GMeanFlowSource(GbyNuMult);
1945  GbyNuMult.clear();
1946 
1947  // Contributions from divU
1948  tmp<volScalarField> divUMult
1949  (
1950  (2.0/3.0)*(zeroFirstCell_*wa()*omega()*gamma_ + ka()*k())
1951  );
1952  meanFlowSource += divUMeanFlowSource(divUMult);
1953 
1954  // Contributions from S2, existing in nut
1955  const volVectorField& U = primalVars_.U();
1956  const volVectorField& Ua = adjointVars_.UaInst();
1957  tmp<volScalarField> nutMeanFlowSourceMult
1958  (
1959  // nut in the diffusion coefficients
1962  + dNutdbMult(U, Ua, nutRef(), "coeffsDiff")
1963  // nut in G
1965  );
1966  meanFlowSource +=
1968  (
1969  nutMeanFlowSourceMult,
1970  F2_,
1971  S_,
1972  case_1_nut_,
1973  gradU_
1974  );
1975 
1976  // G at the first cell includes mag(U.snGrad())
1977  // Add term here
1978  forAll(omega().boundaryFieldRef(), patchi)
1979  {
1980  fvPatchScalarField& omegaWall = omega().boundaryFieldRef()[patchi];
1981  if (isA<omegaWallFunctionFvPatchScalarField>(omegaWall))
1982  {
1983  const fvPatch& patch = mesh_.boundary()[patchi];
1984 
1985  const autoPtr<incompressible::turbulenceModel>& turbModel =
1987 
1988  const scalarField& y = turbModel().y()[patchi];
1989  const tmp<scalarField> tnuw = turbModel().nu(patchi);
1990  const scalarField& nuw = tnuw();
1991 
1992  const nutWallFunctionFvPatchScalarField& nutw =
1993  refCast<nutWallFunctionFvPatchScalarField>
1994  (nutRef().boundaryFieldRef()[patchi]);
1995  const wallFunctionCoefficients& wallCoeffs = nutw.wallCoeffs();
1996  const scalar Cmu = wallCoeffs.Cmu();
1997  const scalar kappa = wallCoeffs.kappa();
1998  const scalar Cmu25 = pow025(Cmu);
1999 
2000  const labelList& faceCells = patch.faceCells();
2001 
2002  const fvPatchVectorField& Uw =
2003  primalVars_.U().boundaryField()[patchi];
2004  vectorField snGradUw(Uw.snGrad());
2005  const scalarField& deltaCoeffs = patch.deltaCoeffs();
2006  forAll(faceCells, facei)
2007  {
2008  const label celli = faceCells[facei];
2009  // Volume will be added when meanFlowSource is added to UaEqn
2010  meanFlowSource[celli] +=
2011  ka()[celli]*case_1_Pk_[celli]
2012  *(nutw[facei] + nuw[facei])
2013  *snGradUw[facei].normalise()
2014  *Cmu25*sqrt(k()[celli])
2015  *deltaCoeffs[facei]
2016  /(kappa*y[facei]);
2017  }
2018  }
2019  }
2020  return tmeanFlowSource;
2021 }
2022 
2023 
2024 tmp<volScalarField> adjointkOmegaSST::nutJacobianTMVar1() const
2025 {
2026  // Compute dnut_dk anew since the current copy of dnut_dk_
2027  // might not be updated
2028  const volVectorField& U = primalVars_.U();
2029  tmp<volScalarField> S2
2030  (
2031  2*magSqr(symm(fvc::grad(U)))
2033  );
2034  volScalarField S(sqrt(S2));
2035  volScalarField F2(this->F2());
2036 
2037  // Computation of nut switches
2038  volScalarField nut_C1(a1_*omega() - b1_*F2*S);
2039  volScalarField arg2_C2
2040  (
2041  (scalar(2)/betaStar_)*sqrt(k())/(omega()*y_)
2042  - scalar(500)*nu()/(sqr(y_)*omega())
2043  );
2044  volScalarField arg2_C1
2045  (
2046  max
2047  (
2048  (scalar(2)/betaStar_)*sqrt(k())/(omega()*y_),
2049  scalar(500)*nu()/(sqr(y_)*omega())
2050  ) - scalar(100)
2051  );
2052  volScalarField case_2_nut(neg0(nut_C1)*pos(arg2_C2)*neg(arg2_C1));
2053 
2054  return dnut_dk(F2, S, case_2_nut);
2055 }
2056 
2057 
2059 {
2060  // Compute dnut_omega anew since the current copy of dnut_domega_
2061  // might not be updated
2062  const volVectorField& U = primalVars_.U();
2064  (
2065  2*magSqr(symm(fvc::grad(U)))
2067  );
2068  volScalarField S(sqrt(S2));
2069  volScalarField F2(this->F2());
2070 
2071  // Computation of nut switches
2072  volScalarField nut_C1(a1_*omega() - b1_*F2*S);
2073  volScalarField arg2_C2
2074  (
2075  (scalar(2)/betaStar_)*sqrt(k())/(omega()*y_)
2076  - scalar(500)*nu()/(sqr(y_)*omega())
2077  );
2078  volScalarField arg2_C1
2079  (
2080  max
2081  (
2082  (scalar(2)/betaStar_)*sqrt(k())/(omega()*y_),
2083  scalar(500)*nu()/(sqr(y_)*omega())
2084  ) - scalar(100)
2085  );
2086  volScalarField case_1_nut(pos(nut_C1));
2087  volScalarField case_2_nut(neg0(nut_C1)*pos(arg2_C2)*neg(arg2_C1));
2088  volScalarField case_3_nut(neg0(nut_C1)*neg0(arg2_C2)*neg(arg2_C1));
2089 
2090  return dnut_domega(F2, S, case_1_nut, case_2_nut, case_3_nut);
2091 }
2092 
2093 
2095 (
2096  tmp<volScalarField>& dNutdUMult
2097 ) const
2098 {
2099  const volVectorField& U = primalVars_.U();
2100  volTensorField gradU(fvc::grad(U));
2102  (
2103  2*magSqr(symm(gradU))
2105  );
2106  volScalarField S(sqrt(S2));
2107  volScalarField F2(this->F2());
2108 
2109  // Computation of nut switches
2110  volScalarField nut_C1(a1_*omega() - b1_*F2*S);
2111  volScalarField case_1_nut(pos(nut_C1));
2112 
2113  return nutMeanFlowSource(dNutdUMult, F2, S, case_1_nut, gradU);
2114 }
2115 
2116 
2118 {
2119  return
2120  (
2121  alphaK_.boundaryField()[patchI]*nutRef().boundaryField()[patchI]
2122  + nu()().boundaryField()[patchI]
2123  );
2124 }
2125 
2126 
2128 {
2129  return
2130  (
2131  alphaOmega_.boundaryField()[patchI]*nutRef().boundaryField()[patchI]
2132  + nu()().boundaryField()[patchI]
2133  );
2134 }
2135 
2136 
2138 {
2140 
2141  if (!adjointTurbulence_)
2142  {
2143  return;
2144  }
2145 
2147 
2148  // Primal and adjoint fields
2149  const volVectorField& U = primalVars_.U();
2151  volScalarField dR_dnut(this->dR_dnut());
2153 
2155 
2156  tmp<fvScalarMatrix> waEqn
2157  (
2158  fvm::div(-phi, wa())
2162  + waEqnSourceFromF1()
2164  + fvm::Sp(zeroFirstCell_*scalar(2.)*beta_*omega(), wa())
2165  + fvm::SuSp
2166  (
2167  zeroFirstCell_()*gamma_*((2.0/3.0)*divU - dGPrime_domega().ref()()),
2168  wa()
2169  )
2170  - (case_2_Pk_*c1_ - scalar(1))*betaStar_*k()*ka()
2171  ==
2172  fvOptions(wa())
2173  );
2174 
2175  // Boundary manipulate changes the diagonal component, so relax has to
2176  // come after that
2177  waEqn.ref().boundaryManipulate(wa().boundaryFieldRef());
2178  waEqn.ref().relax();
2179 
2180  // Sources from the objective should be added after the boundary
2181  // manipulation
2182  objectiveManager_.addSource(waEqn.ref());
2183  fvOptions.constrain(waEqn.ref());
2184  waEqn.ref().solve();
2185 
2186  // Adjoint Turbulent kinetic energy equation
2187  tmp<fvScalarMatrix> kaEqn
2188  (
2189  fvm::div(-phi, ka())
2190  + fvm::SuSp(fvc::div(phi), ka())
2191  - fvm::laplacian(DkEff_, ka())
2192  + fvm::Sp(betaStar_*omega(), ka())
2193  - case_2_Pk_()*c1_*betaStar_*omega()()*ka()()
2194  + fvm::SuSp(scalar(2.0/3.0)*divU, ka())
2196  + kaEqnSourceFromF1()
2197  + dR_dnut()*dnut_dk_()
2198  - zeroFirstCell_()*gamma_*dGPrime_dk().ref()()*wa()
2199  ==
2200  fvOptions(ka())
2201  );
2202 
2203  kaEqn.ref().relax();
2204  kaEqn.ref().boundaryManipulate(ka().boundaryFieldRef());
2205  addWallFunctionTerms(kaEqn.ref(), dR_dnut);
2206  fvOptions.constrain(kaEqn.ref());
2207  // Add sources from the objective functions
2208  objectiveManager_.addSource(kaEqn.ref());
2209 
2210  kaEqn.ref().solve();
2211 
2213  {
2214  dimensionedScalar maxwa = max(mag(wa()));
2216  Info<< "Max mag (" << wa().name() << ") = " << maxwa.value() << endl;
2217  Info<< "Max mag (" << ka().name() << ") = " << maxka.value() << endl;
2218  }
2219 }
2220 
2223 {
2224  return adjMomentumBCSourcePtr_();
2225 }
2226 
2227 
2229 {
2232 
2233  forAll(mesh_.boundary(), patchi)
2234  {
2235  vectorField nf(mesh_.boundary()[patchi].nf());
2236  wallShapeSens[patchi] = nf & FITerm.boundaryField()[patchi];
2237  }
2238  return wallShapeSens;
2239 }
2240 
2243 {
2244  return wallFloCoSensitivitiesPtr_();
2245 }
2246 
2247 
2249 {
2251  (
2252  IOobject
2253  (
2254  "adjointEikonalSource" + type(),
2255  runTime_.timeName(),
2256  mesh_,
2259  ),
2260  mesh_,
2262  );
2263 }
2264 
2265 
2267 {
2268  const volVectorField& U = primalVars_.U();
2269  const volScalarField& kInst =
2270  primalVars_.RASModelVariables()->TMVar1Inst();
2271  const volScalarField& omegaInst =
2272  primalVars_.RASModelVariables()->TMVar2Inst();
2273 
2274  tmp<volScalarField> arg1 = min
2275  (
2276  min
2277  (
2278  max
2279  (
2280  (scalar(1)/betaStar_)*sqrt(k())/(omega()*y_),
2281  scalar(500)*nu()/(sqr(y_)*omega())
2282  ),
2283  (4*alphaOmega2_)*k()/(CDkOmegaPlus_*sqr(y_))
2284  ),
2285  scalar(10)
2286  );
2287 
2288  // Interpolation schemes used by the primal convection terms
2289  auto kScheme(convectionScheme(kInst.name()));
2290  auto omegaScheme(convectionScheme(omegaInst.name()));
2291  const surfaceVectorField& Sf = mesh_.Sf();
2292  tmp<volTensorField> tFISens
2293  (
2295  (
2296  IOobject
2297  (
2298  type() + "FISensTerm",
2299  mesh_.time().timeName(),
2300  mesh_,
2303  ),
2304  mesh_,
2307  )
2308  );
2309  volTensorField& FISens = tFISens.ref();
2310  FISens =
2311  // k convection
2312  - ka()*fvc::div
2313  (
2314  kScheme().interpolate(k())
2316  )
2317  // k diffusion
2318  + ka()*T(fvc::grad(DkEff_*gradK_))
2319  - DkEff_*(fvc::grad(ka())*gradK_)
2320  // omega convection
2322  (
2323  omegaScheme().interpolate(omega())
2325  )
2326  // omega diffusion
2329  // terms including GbyNu0
2330  + (
2332  + case_1_Pk_*ka()*nutRef()
2334  // S2 (includes contribution from nut in UEqn as well)
2335  + (
2336  dR_dnut()*a1_*k()/(b1_*S_*S_*S_*F2_)
2338  *(c1_/a1_)*betaStar_*omega()*b1_*F2_/S_
2339  )*T(gradU_ & twoSymm(gradU_))*(1. - case_1_nut_)
2340  // CDkOmega in omegaEqn
2341  + 2.*wa()*(1. - F1_)*alphaOmega2_/omega()*zeroFirstCell_
2343  // F1
2344  - dR_dF1()
2345  *(dF1_dGradK(arg1)*gradK_ + dF1_dGradOmega(arg1)*gradOmega_);
2346 
2348 
2349  return tFISens;
2350 }
2351 
2352 
2354 (
2355  const word& designVarsName
2356 ) const
2357 {
2359  auto tres(tmp<scalarField>::New(mesh_.nCells(), Zero));
2360 
2361  // Sensitivity from the source term in the k equation
2362  scalarField auxSens
2363  (k().primitiveField()*ka().primitiveField());
2365  (
2366  tres.ref(), auxSens, fvOptions, k().name(), designVarsName
2367  );
2368 
2369  // Sensitivity from the source term in the omega equation
2370  auxSens = omega().primitiveField()*wa().primitiveField();
2372  (
2373  tres.ref(), auxSens, fvOptions, omega().name(), designVarsName
2374  );
2375 
2376  return tres;
2377 }
2378 
2379 
2381 {
2384 }
2385 
2386 
2388 {
2389  if (adjointRASModel::read())
2390  {
2392  alphaK1_.readIfPresent(this->coeffDict());
2393  alphaK2_.readIfPresent(this->coeffDict());
2396  gamma1_.readIfPresent(this->coeffDict());
2397  gamma2_.readIfPresent(this->coeffDict());
2398  beta1_.readIfPresent(this->coeffDict());
2399  beta2_.readIfPresent(this->coeffDict());
2401  a1_.readIfPresent(this->coeffDict());
2402  b1_.readIfPresent(this->coeffDict());
2403  c1_.readIfPresent(this->coeffDict());
2404  F3_.readIfPresent("F3", this->coeffDict());
2405 
2406  return true;
2407  }
2408  else
2409  {
2410  return false;
2411  }
2412 }
2413 
2414 
2415 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
2416 
2417 } // End namespace adjointRASModels
2418 } // End namespace incompressible
2419 } // End namespace Foam
2420 
2421 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
volScalarField DOmegaEff_
Diffusivity of the omega equation.
virtual tmp< volScalarField > GbyNu(const volScalarField &GbyNu0, const volScalarField &F2, const volScalarField &S2) const
tmp< volVectorField > nutMeanFlowSource(tmp< volScalarField > &mult, const volScalarField &F2, const volScalarField &S, const volScalarField &case_1_nut, const volTensorField &gradU) const
Contributions from nut(U)
dimensionedScalar tanh(const dimensionedScalar &ds)
virtual const boundaryVectorField & wallShapeSensitivities()
Sensitivity derivative contributions when using the (E)SI approach.
const Type & value() const noexcept
Return const reference to value.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
bool includeDistance_
Does the turbulence model include distances and should the adjoint to the distance field be computed...
tmp< volScalarField > dnut_dk(const volScalarField &F2, const volScalarField &S, const volScalarField &case_2_nut) const
Nut Jacobian wrt k.
virtual tmp< scalarField > topologySensitivities(const word &designVarsName) const
Term contributing to the computation of topology optimisation sensitivities.
const surfaceVectorField & Sf() const
Return cell face area vectors.
tmp< volVectorField > dF1_dGradK(const volScalarField &arg1) const
F1 Jacobian wrt grad(k)
fvPatchField< vector > fvPatchVectorField
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:213
dimensionedScalar log(const dimensionedScalar &ds)
Inversed weight central-differencing interpolation scheme class.
Definition: reverseLinear.H:53
tmp< volScalarField > dF2_domega(const volScalarField &F2, const volScalarField &case_2_nut, const volScalarField &case_3_nut) const
F2 Jacobian wrt omega.
tmp< volScalarField > diffusionNutMeanFlowMult(const volScalarField &primalField, const volScalarField &adjointField, const volScalarField &coeffField) const
Contributions from nut(U), in the diffusion coefficients of the turbulence model. ...
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:160
tmp< volVectorField > dF1_dGradOmega(const volScalarField &arg1) const
F1 Jacobian wrt grad(omega)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
tmp< volScalarField > waEqnSourceFromF1() const
Source to waEqn from the differentiation of F1.
Central-differencing interpolation scheme class.
Definition: linear.H:51
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
const dimensionedScalar G
Newtonian constant of gravitation.
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
const wallFunctionCoefficients & wallCoeffs() const noexcept
Return wallFunctionCoefficients.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
bool readIfPresent(const dictionary &dict)
Update the value of dimensioned<Type> if found in the dictionary, lookup in dictionary with the name(...
Class for managing objective functions.
virtual bool read()
Read adjointRASProperties dictionary.
virtual tmp< volSymmTensorField > devReff() const
Return the effective stress tensor including the laminar stress.
volScalarField dnut_domega_
Nut Jacobian w.r.t. omega.
static tmp< edgeInterpolationScheme< Type > > scheme(const edgeScalarField &faceFlux, Istream &schemeData)
Return weighting factors for scheme given from Istream.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
tmp< volScalarField > kaEqnSourceFromF1() const
Source to kaEqn from the differentiation of F1.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
dimensionedScalar pow025(const dimensionedScalar &ds)
virtual const boundaryVectorField & wallFloCoSensitivities()
Sensitivity terms for flow control, emerging from the turbulence model differentiation.
tmp< volScalarField > dF1_domega(const volScalarField &arg1) const
F1 Jacobian wrt omega (no contributions from grad(omega))
tmp< volVectorField > GMeanFlowSource(tmp< volSymmTensorField > &GbyNuMult) const
Contributions from the G.
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, any/none. Also accepts 0/1 as a string and shortcuts t/f, y/n.
Definition: Switch.H:77
static void nullifyField(GeometricField< Type, PatchField, GeoMesh > &fieldPtr)
Nullify field and old times, if present.
#define F2(B, C, D)
Definition: SHA1.C:150
const surfaceScalarField & phiInst() const
Return const reference to volume flux.
autoPtr< boundaryVectorField > adjMomentumBCSourcePtr_
Source to the adjoint momentum BC emerging from differentiating the turbulence model.
tmp< volScalarField::Internal > beta(const volScalarField::Internal &F1) const
virtual tmp< volVectorField > adjointMeanFlowSource()
Source term added to the adjoint mean flow due to the differentiation of the turbulence model...
Generic dimensioned Type class.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
autoPtr< boundaryVectorField > wallShapeSensitivitiesPtr_
Wall sensitivity term for shape optimisation.
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:52
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
tmp< volScalarField > kaEqnSourceFromCDkOmega() const
Source to kaEqn from the differentiation of CDkOmega.
Finite-volume options.
Definition: fvOptions.H:51
virtual tmp< volScalarField > nutJacobianTMVar1() const
Jacobian of nut wrt to k.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
virtual tmp< volScalarField > distanceSensitivities()
Contributions to the adjoint eikonal equation (zero for now)
tmp< volScalarField > alphaK(const volScalarField &F1) const
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
tmp< volScalarField > alphaOmega(const volScalarField &F1) const
dimensionedScalar neg(const dimensionedScalar &ds)
tmp< volVectorField > convectionMeanFlowSource(const volScalarField &primalField, const volScalarField &adjointField) const
Contributions from the turbulence model convection terms.
Macros for easy insertion into run-time selection tables.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
fv::options & fvOptions
tmp< volScalarField > nu() const
Return the laminar viscosity.
tmp< volScalarField > DkEff(const volScalarField &F1) const
Return the effective diffusivity for k.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
void addWallFunctionTerms(fvScalarMatrix &kaEqn, const volScalarField &dR_dnut)
Contributions from the differentiation of k existing in nutkWallFunction.
Switch adjointTurbulence_
Turbulence on/off flag.
tmp< surfaceInterpolationScheme< scalar > > convectionScheme(const word &varName) const
Return the interpolation scheme used by the primal convection term of the equation corresponding to t...
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
autoPtr< boundaryVectorField > wallFloCoSensitivitiesPtr_
Wall sensitivity term for flow control optimisation.
word timeName
Definition: getTimeIndex.H:3
dimensionedScalar pos(const dimensionedScalar &ds)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
const surfaceScalarField & phi() const
Return const reference to volume flux.
Base class for solution control classes.
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
bool changedPrimalSolution_
Has the primal solution changed?
The class nutWallFunction is an abstract base class that hosts calculation methods and common functi...
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
void setSize(const label n)
Alias for resize()
Definition: List.H:316
virtual tmp< volScalarField > nutJacobianTMVar2() const
Jacobian of nut wrt to omega.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
virtual void correct()
Solve the adjoint turbulence equations.
bool useSolverNameForFields() const
Append solver name to fields?
Definition: variablesSet.C:83
virtual tmp< volTensorField > FISensitivityTerm()
Sensitivity derivative contributions when using the FI approach.
tmp< volScalarField > DomegaEff(const volScalarField &F1) const
Return the effective diffusivity for omega.
fvPatchField< scalar > fvPatchScalarField
volTensorField gradU_
Cached primal gradient fields.
A class for handling words, derived from Foam::string.
Definition: word.H:63
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
virtual tmp< fvVectorMatrix > divDevReff(volVectorField &U) const
Return the transpose part of the adjoint momentum stresses.
Abstract base class for incompressible turbulence models.
bool printMaxMags() const
Print max mags of solver fields.
dimensionedScalar neg0(const dimensionedScalar &ds)
GeometricField< tensor, fvsPatchField, surfaceMesh > surfaceTensorField
tmp< volVectorField > divUMeanFlowSource(tmp< volScalarField > &divUMult) const
Contributions from the divU.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
tmp< volScalarField > dR_dnut()
Derivative of the primal equations wrt nut.
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
wordList adjointTMVariablesBaseNames_
Base names of the adjoint fields.
addToRunTimeSelectionTable(adjointRASModel, adjointkOmegaSST, dictionary)
Calculation of adjoint based sensitivities for topology optimisation. This returns just the field par...
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
autoPtr< volScalarField > adjointTMVariable1Ptr_
Adjoint turbulence model variable 1.
virtual void nullify()
Nullify all adjoint turbulence model fields and their old times.
const word & solverName() const
Return solver name.
Definition: variablesSet.C:77
virtual void correct()
Solve the adjoint turbulence equations.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
virtual const tmp< volScalarField > nut() const
Return the turbulence viscosity.
autoPtr< volScalarField > adjointTMVariable2Ptr_
Adjoint turbulence model variable 2.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
tmp< volScalarField > dNutdbMult(const volScalarField &primalField, const volScalarField &adjointField, const volScalarField &coeffField, const volScalarField &bcField, const word &schemeName) const
Term multiplying dnut/db, coming from the turbulence model.
virtual tmp< volVectorField > nutJacobianU(tmp< volScalarField > &dNutdUMult) const
Jacobian of nut wrt the flow velocity.
Us boundaryFieldRef().evaluateCoupled< coupledFaPatch >()
const volVectorField & U() const
Return const reference to velocity.
volScalarField DkEff_
Diffusivity of the k equation.
volScalarField case_1_Pk_
Switch fields for the production in the k Eqn.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
dimensionedScalar pos0(const dimensionedScalar &ds)
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
Definition: fvPatchField.C:221
rDeltaT ref()
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
tmp< volScalarField > dF1_dk(const volScalarField &arg1) const
F1 Jacobian wrt k (no contributions from grad(k))
UPtrList< const labelUList > faceCells() const
Return a list of faceCells for each patch.
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1101
scalar Cmu() const noexcept
Return the object: Cmu.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const solverControl & getSolverControl() const
Return const reference to solverControl.
tmp< volScalarField > dF2_dk(const volScalarField &F2, const volScalarField &case_2_nut) const
F2 Jacobian wrt k.
static void postProcessSens(scalarField &sens, scalarField &auxSens, fv::options &fvOptions, const word &fieldName, const word &designVariablesName)
Add part of the sensitivities coming from fvOptions.
virtual void addSource(fvVectorMatrix &matrix)
Add contribution to adjoint momentum PDEs.
tmp< fvScalarMatrix > waEqnSourceFromCDkOmega() const
Source to waEqn from the differentiation of CDkOmega.
const autoPtr< incompressible::RASModelVariables > & RASModelVariables() const
Return const reference to the turbulence model variables.
Manages the adjoint mean flow fields and their mean values.
Continuous adjoint to the kOmegaSST turbulence model for incompressible flows.
const volVectorField & UaInst() const
Return const reference to velocity.
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
Field< Type > & source() noexcept
Definition: fvMatrix.H:533
virtual tmp< scalarField > diffusionCoeffVar1(label patchI) const
Diffusion coeff at the boundary for k.
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:183
U
Definition: pEqn.H:72
tmp< volScalarField > coeffsDifferentiation(const volScalarField &primalField, const volScalarField &adjointField, const word &schemeName) const
Differentiation of the turbulence model diffusion coefficients.
This boundary condition provides a wall function for the specific dissipation rate (i...
ITstream & divScheme(const word &name) const
Get div scheme for given name, or default.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
void updatePrimalRelatedFields()
Update of the primal cached fields.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label nCells() const noexcept
Number of mesh cells.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
scalar yPlus
IOobject(const IOobject &)=default
Copy construct.
tmp< volScalarField > dnut_domega(const volScalarField &F2, const volScalarField &S, const volScalarField &case_1_nut, const volScalarField &case_2_nut, const volScalarField &case_3_nut) const
Nut Jacobian wrt omega.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
dimensionedScalar pow4(const dimensionedScalar &ds)
Nothing to be read.
const std::string patch
OpenFOAM patch number as a std::string.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
tmp< volScalarField::Internal > gamma(const volScalarField::Internal &F1) const
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
virtual tmp< scalarField > diffusionCoeffVar2(label patchI) const
Diffusion coeff at the boundary for omega.
messageStream Info
Information stream (stdout output on master, null elsewhere)
objectiveManager & objectiveManager_
Reference to the objectiveManager.
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
const autoPtr< incompressible::turbulenceModel > & turbulence() const
Return const reference to the turbulence model.
Field< vector > vectorField
Specialisation of Field<T> for vector.
virtual tmp< volVectorField > nonConservativeMomentumSource() const
Non-conservative part of the terms added to the mean flow equations.
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition: tmpI.H:289
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
void boundaryManipulate(typename GeometricField< Type, fvPatchField, volMesh >::Boundary &values)
Manipulate based on a boundary field.
Definition: fvMatrix.C:1267
scalar kb
tmp< volScalarField > dGPrime_dk() const
GbyNu Jacobian wrt k.
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
zeroField divU
Definition: alphaSuSp.H:3
List< label > labelList
A List of labels.
Definition: List.H:62
A class for managing temporary objects.
Definition: HashPtrTable.H:50
volScalarField S2_
Primal cached fields involved in the solution of the.
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:96
Class to host the wall-function coefficients being used in the wall function boundary conditions...
virtual bool read()
Read adjointRASProperties dictionary.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
const nearWallDist & y() const
Return the near wall distances.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
tmp< volScalarField > dGPrime_domega() const
GbyNu Jacobian wrt omega.
virtual const boundaryVectorField & adjointMomentumBCSource() const
Source for the outlet adjoint momentum BC coming from differentiating the turbulence model...
bool readIfPresent(const word &key, const dictionary &dict)
Update the value of the Switch if it is found in the dictionary.
Definition: Switch.C:342
#define M(I)
tmp< volScalarField > dR_dF1() const
Derivative of the primal equations wrt F1.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:40
An input stream of tokens.
Definition: ITstream.H:52
scalar nut
static void setField(autoPtr< GeometricField< Type, fvPatchField, volMesh >> &fieldPtr, const fvMesh &mesh, const word &baseName, const word &solverName, const bool useSolverNameForFields)
Read vol fields.
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:56
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
Definition: SymmTensorI.H:491
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const dimensionSet dimVelocity