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 = volScalarField::New
129  (
130  "zeroFirstCell",
132  mesh_,
134  );
135  auto& zeroFirstCell = tzeroFirstCell.ref();
136 
138  label counter(0);
139  for (const fvPatchScalarField& omegab : omega().boundaryField())
140  {
141  const fvPatch& patch = omegab.patch();
142  if (isA<omegaWallFunctionFvPatchScalarField>(omegab))
143  {
144  const label patchi = patch.index();
145  const labelList& faceCells = patch.faceCells();
146  fvPatchScalarField& bf = zeroFirstCell.boundaryFieldRef()[patchi];
147  forAll(faceCells, faceI)
148  {
149  const label cellI = faceCells[faceI];
150 
151  zeroFirstCell[cellI] = 0.;
152  bf[faceI] = 0.;
153  firstCellIDs_[counter++] = cellI;
154  }
155  }
156  }
157  firstCellIDs_.setSize(counter);
159  zeroFirstCell.correctBoundaryConditions();
160 
161  return tzeroFirstCell;
162 }
163 
164 
166 {
167  const volVectorField& U = primalVars_.U();
168  const volVectorField& Ua = adjointVars_.UaInst();
169  word scheme("coeffsDiff");
170  tmp<volScalarField> tdR_dnut
171  (
172  dNutdbMult(U, Ua, nutRef(), scheme)
173  + dNutdbMult(k(), ka(), alphaK_, nutRef(), scheme)
176  );
177  volScalarField& dRdnut = tdR_dnut.ref();
178 
179  forAll(mesh_.boundary(), pI)
180  {
181  const fvPatch& patch = mesh_.boundary()[pI];
182  const fvPatchScalarField& nutb = nutRef().boundaryField()[pI];
183  const labelList& faceCells = patch.faceCells();
184  if (isA<nutkWallFunctionFvPatchScalarField>(nutb))
185  {
186  fvPatchScalarField& bf = dRdnut.boundaryFieldRef()[pI];
187  const scalarField kSnGrad(k().boundaryField()[pI].snGrad());
188  const scalarField omegaSnGrad
189  (
190  omega().boundaryField()[pI].snGrad()
191  );
192  const fvPatchScalarField& akb = alphaK_.boundaryField()[pI];
193  const fvPatchScalarField& aOmegab = alphaOmega_.boundaryField()[pI];
194  const vectorField USnGrad(U.boundaryField()[pI].snGrad());
195  const fvPatchTensorField& gradUb = gradU_.boundaryField()[pI];
196  const vectorField nf(mesh_.boundary()[pI].nf());
197  forAll(faceCells, fI)
198  {
199  const label cI(faceCells[fI]);
200  bf[fI] =
201  - wa()[cI]*zeroFirstCell_[cI]*aOmegab[fI]*omegaSnGrad[fI]
202  - ka()[cI]*akb[fI]*kSnGrad[fI]
203  - (Ua[cI] & (USnGrad[fI] + (dev2(gradUb[fI]) & nf[fI])));
204  }
205  }
206  }
207 
208  return tdR_dnut;
209 }
210 
211 
213 (
214  const volScalarField& F2,
215  const volScalarField& S,
216  const volScalarField& case_1_nut,
217  const volScalarField& case_2_nut,
218  const volScalarField& case_3_nut
219 ) const
220 {
221  return
222  (
223  - case_1_nut*k()/sqr(omega())
224  - a1_*k()/(b1_*S*F2*F2)*dF2_domega(F2, case_2_nut, case_3_nut)
225  );
226 }
227 
228 
230 (
231  const volScalarField& F2,
232  const volScalarField& S,
233  const volScalarField& case_2_nut
234 ) const
235 {
236  return
237  (
239  - a1_*k()/(b1_*S*F2*F2)*dF2_dk(F2, case_2_nut)
240  );
241 }
242 
243 
245 (
246  const volScalarField& F2,
247  const volScalarField& case_2_nut,
248  const volScalarField& case_3_nut
249 ) const
250 {
251  tmp<volScalarField> arg2 = min
252  (
253  max
254  (
255  (scalar(2)/betaStar_)*sqrt(k())/(omega()*y_),
256  scalar(500)*nu()/(sqr(y_)*omega())
257  ),
258  scalar(100)
259  );
260 
261  return
262  - scalar(2)*arg2*(1 - F2*F2)*
263  (
264  case_2_nut*scalar(2)*sqrt(k())/(betaStar_*sqr(omega())*y_)
265  + case_3_nut*scalar(500)*nu()/sqr(omega()*y_)
266  );
267 }
268 
269 
271 (
272  const volScalarField& F2,
273  const volScalarField& case_2_nut
274 ) const
275 {
276  tmp<volScalarField> arg2 = min
277  (
278  max
279  (
280  (scalar(2)/betaStar_)*sqrt(k())/(omega()*y_),
281  scalar(500)*nu()/(sqr(y_)*omega())
282  ),
283  scalar(100)
284  );
286  return
287  case_2_nut*scalar(2)*arg2*(1 - F2*F2)
288  /(betaStar_*omega()*y_*sqrt(k()));
289 }
290 
291 
293 {
294  return
295  (
297  *(
298  max(a1_*omega(), b1_*F2_*S_)
299  + case_1_nut_*a1_*omega()
300  + (scalar(1) - case_1_nut_)*omega()*b1_*S_
302  )
303  );
304 }
305 
306 
308 {
309  return
312 }
313 
314 
316 {
317  const volVectorField& U = primalVars_.U();
319  tmp<volScalarField> tGPrime = GbyNu(GbyNu0_, F2_, S2_);
321  word scheme("coeffsDiff");
322 
323  tmp<volScalarField> tdRdF1
324  (
325  nutRef()
326  *(
330  )
331  );
332  volScalarField& dRdF1 = tdRdF1.ref();
333 
334  typedef fixedValueFvPatchScalarField fixedValue;
335  typedef zeroGradientFvPatchScalarField zeroGrad;
336  typedef omegaWallFunctionFvPatchScalarField omegaWF;
337  forAll(mesh_.boundary(), pI)
338  {
339  const fvPatchScalarField& kb = k().boundaryField()[pI];
340  const fvPatchScalarField& omegab = omega().boundaryField()[pI];
341  fvPatchScalarField& dRdF1b = dRdF1.boundaryFieldRef()[pI];
342  if
343  (
344  isA<zeroGrad>(kb)
345  && (isA<zeroGrad>(omegab) || isA<omegaWF>(omegab))
346  )
347  {
348  dRdF1b = dRdF1b.patchInternalField();
349  }
350  else if (isA<fixedValue>(kb) && isA<fixedValue>(omegab)
351  )
352  {
353  // Note: might need revisiting
354  dRdF1b = dRdF1b.patchInternalField();
355  }
356  }
357 
358  volScalarField dR_dF1_coeffs
359  (
361  *(
362  (gamma1_ - gamma2_)*((2.0/3.0)*omega()*tdivU - tGPrime)
363  + omega()*omega()*(beta1_ - beta2_) + CDkOmega_
364  )
365  );
366 
367  forAll(mesh_.boundary(), pI)
368  {
369  const fvPatchScalarField& kb = k().boundaryField()[pI];
370  const fvPatchScalarField& omegab = omega().boundaryField()[pI];
371  fvPatchScalarField& dRdF1b = dR_dF1_coeffs.boundaryFieldRef()[pI];
372  if
373  (
374  isA<zeroGrad>(kb)
375  && (isA<zeroGrad>(omegab) || isA<omegaWF>(omegab))
376  )
377  {
378  dRdF1b = Zero;
379  }
380  else if (isA<fixedValue>(kb) && isA<fixedValue>(omegab))
381  {
382  dRdF1b = dRdF1b.patchInternalField();
383  }
384  }
386  dRdF1 += dR_dF1_coeffs;
387  return tdRdF1;
388 }
389 
390 
392 (
393  const volScalarField& arg1
394 ) const
395 {
396  return
397  (
398  4*pow3(arg1)*(scalar(1) - F1_*F1_)
399  *(
400  - case_1_F1_*sqrt(k())/(betaStar_*omega()*omega()*y_)
401  - case_2_F1_*scalar(500)*nu()/sqr(omega()*y_)
404  )
405  );
406 }
407 
408 
410 (
411  const volScalarField& arg1
412 ) const
413 {
414  return
415  (
416  - case_3_F1_*scalar(4)*pow3(arg1)*(scalar(1) - F1_*F1_)
417  *scalar(8)*k()*sqr(alphaOmega2_/(CDkOmegaPlus_*y_))/omega()*gradK_
418  );
419 }
420 
421 
423 {
424  tmp<volScalarField> arg1 = min
425  (
426  min
427  (
428  max
429  (
430  (scalar(1)/betaStar_)*sqrt(k())/(omega()*y_),
431  scalar(500)*nu()/(sqr(y_)*omega())
432  ),
434  ),
435  scalar(10)
436  );
437 
438  volScalarField dR_dF1(this->dR_dF1());
439  volScalarField dF1_domega(this->dF1_domega(arg1));
441 
442  surfaceScalarField dR_dGradOmega
443  (
444  interpolationScheme<vector>("div(dR_dGradOmega)")().
446  );
447 
448  return
449  (
451  - fvc::div(dR_dGradOmega)
452  );
453 }
454 
455 
457 {
458  tmp<volVectorField> tsource
459  (
461  );
462  volVectorField& source = tsource.ref();
463 
464  forAll(mesh_.boundary(), pI)
465  {
466  const fvPatchScalarField& omegab = omega().boundaryField()[pI];
467  fvPatchVectorField& sourceb = source.boundaryFieldRef()[pI];
468  if
469  (
470  isA<zeroGradientFvPatchScalarField>(omegab)
471  || isA<omegaWallFunctionFvPatchScalarField>(omegab)
472  )
473  {
474  sourceb = Zero;
475  }
476  else if (isA<fixedValueFvPatchScalarField>(omegab))
477  {
478  sourceb = sourceb.patchInternalField();
479  }
480  }
481 
482  surfaceScalarField sourceFaces
483  (
484  interpolationScheme<vector>("sourceAdjontkOmegaSST")().
485  interpolate(source) & mesh_.Sf()
486  );
487 
488  return
489  (
490  // Differentiation of omega in CDkOmega
492  // Differentiation of grad(omega) in CDkOmega
493  + fvc::div(sourceFaces)
494  );
495 }
496 
497 
499 {
500  tmp<volVectorField> tsource
501  (
503  );
504  volVectorField& source = tsource.ref();
505 
506  forAll(mesh_.boundary(), pI)
507  {
508  const fvPatchScalarField& kb = k().boundaryField()[pI];
509  fvPatchVectorField& sourceb = source.boundaryFieldRef()[pI];
510  if (isA<zeroGradientFvPatchScalarField>(kb))
511  {
512  sourceb = Zero;
513  }
514  else if (isA<fixedValueFvPatchScalarField>(kb))
515  {
516  sourceb = sourceb.patchInternalField();
517  }
518  }
519 
520  return
521  fvc::div
522  (
523  interpolationScheme<vector>("sourceAdjontkOmegaSST")().
524  interpolate(source) & mesh_.Sf()
525  );
526 }
527 
528 
530 (
531  const volScalarField& arg1
532 ) const
533 {
534  return
535  (
536  4*pow3(arg1)*(scalar(1) - F1_*F1_)
537  *(
538  case_1_F1_*0.5/(betaStar_*omega()*sqrt(k())*y_)
540  )
541  );
542 }
543 
544 
546 (
547  const volScalarField& arg1
548 ) const
549 {
550  return
551  (
552  - case_3_F1_*scalar(4)*pow3(arg1)*(scalar(1) - F1_*F1_)
553  *scalar(8)*k()*sqr(alphaOmega2_/(CDkOmegaPlus_*y_))/omega()
554  *gradOmega_
555  );
556 }
557 
558 
560 {
561  tmp<volScalarField> arg1 = min
562  (
563  min
564  (
565  max
566  (
567  (scalar(1)/betaStar_)*sqrt(k())/(omega()*y_),
568  scalar(500)*nu()/(sqr(y_)*omega())
569  ),
571  ),
572  scalar(10)
573  );
574 
575  volScalarField dR_dF1(this->dR_dF1());
576  volScalarField dF1_dk(this->dF1_dk(arg1));
577  volVectorField dF1_dGradK(this->dF1_dGradK(arg1));
578 
579  surfaceScalarField dR_dGradK
580  (
581  interpolationScheme<vector>("div(dR_dGradK)")().
583  );
584 
585  return
586  (
588  - fvc::div(dR_dGradK)
589  );
590 }
591 
592 
594 (
595  const volScalarField& primalField,
596  const volScalarField& adjointField,
597  const word& schemeName
598 ) const
599 {
601  (interpolationScheme<scalar>(schemeName));
602  surfaceScalarField snGradPrimal(fvc::snGrad(primalField)*mesh_.magSf());
603  surfaceScalarField flux(scheme().interpolate(adjointField)*snGradPrimal);
604 
605  forAll(mesh_.boundary(), pI)
606  {
607  const fvPatchScalarField& primalbf = primalField.boundaryField()[pI];
608  if (isA<zeroGradientFvPatchScalarField>(primalbf))
609  {
610  // Needless, but just to be safe
611  snGradPrimal.boundaryFieldRef()[pI] = Zero;
612  flux.boundaryFieldRef()[pI] = Zero;
613  }
614  else if (isA<fixedValueFvPatchScalarField>(primalbf))
615  {
616  // Note: to be potentially revisited
617  snGradPrimal.boundaryFieldRef()[pI] = Zero;
618  flux.boundaryFieldRef()[pI] = Zero;
619  }
620  }
621 
622  return (fvc::div(flux) - adjointField*fvc::div(snGradPrimal));
623 }
624 
625 
627 (
628  const volScalarField& primalField,
629  const volScalarField& adjointField,
630  const volScalarField& coeffField,
631  const volScalarField& bcField,
632  const word& schemeName
633 ) const
634 {
636  (interpolationScheme<scalar>(schemeName));
637  surfaceScalarField snGradPrimal(fvc::snGrad(primalField)*mesh_.magSf());
638  surfaceScalarField flux(scheme().interpolate(adjointField)*snGradPrimal);
639 
640  forAll(mesh_.boundary(), pI)
641  {
642  const fvPatchScalarField& bc = bcField.boundaryField()[pI];
643  if (isA<zeroGradientFvPatchScalarField>(bc))
644  {
645  const fvPatchScalarField& coeffFieldb =
646  coeffField.boundaryField()[pI];
647  snGradPrimal.boundaryFieldRef()[pI] *=
648  coeffFieldb/coeffFieldb.patchInternalField();
649  flux.boundaryFieldRef()[pI] = Zero;
650  }
651  else if (isA<fixedValueFvPatchScalarField>(bc))
652  {
653  snGradPrimal.boundaryFieldRef()[pI] = Zero;
654  flux.boundaryFieldRef()[pI] = Zero;
655  }
656  }
657 
658  return coeffField*(fvc::div(flux) - adjointField*fvc::div(snGradPrimal));
659 }
660 
661 
663 (
664  const volVectorField& U,
665  const volVectorField& Ua,
666  const volScalarField& nut,
667  const word& schemeName
668 ) const
669 {
671  (interpolationScheme<vector>(schemeName));
673  surfaceScalarField flux(scheme().interpolate(Ua) & snGradU);
674 
675  // Terms form the Laplacian part of the momentum stresses
676  forAll(mesh_.boundary(), pI)
677  {
678  const fvPatchScalarField& bc = nut.boundaryField()[pI];
679  if (isA<zeroGradientFvPatchScalarField>(bc))
680  {
681  flux.boundaryFieldRef()[pI] = Zero;
682  }
683  else if (isA<fixedValueFvPatchScalarField>(bc))
684  {
685  snGradU.boundaryFieldRef()[pI] = Zero;
686  flux.boundaryFieldRef()[pI] = Zero;
687  }
688  }
689 
690  // Terms form the tranpose gradient in the momentum stresses
691  const surfaceVectorField& Sf = mesh_.Sf();
692  surfaceTensorField fluxTranspose
693  (
694  reverseLinear<vector>(mesh_).interpolate(Ua)*Sf
695  );
696  forAll(mesh_.boundary(), pI)
697  {
698  const fvPatchVectorField& Ub = U.boundaryField()[pI];
699  if (!isA<coupledFvPatchVectorField>(Ub))
700  {
701  const vectorField Uai(Ua.boundaryField()[pI].patchInternalField());
702  const vectorField& Sfb = Sf.boundaryField()[pI];
703  fluxTranspose.boundaryFieldRef()[pI] = Uai*Sfb;
704  }
705  }
706  volScalarField M(dev2(gradU_) && fvc::div(fluxTranspose));
707  const DimensionedField<scalar, volMesh>& V = mesh_.V();
708  forAll(mesh_.boundary(), pI)
709  {
710  const fvPatchScalarField& bc = nut.boundaryField()[pI];
711  if (isA<zeroGradientFvPatchScalarField>(bc))
712  {
713  const vectorField Uai(Ua.boundaryField()[pI].patchInternalField());
714  const tensorField dev2GradU
715  (
716  dev2(gradU_.boundaryField()[pI].patchInternalField())
717  );
718  const vectorField& Sfb = Sf.boundaryField()[pI];
719  const labelList& faceCells = mesh_.boundary()[pI].faceCells();
720  forAll(faceCells, fI)
721  {
722  const label celli = faceCells[fI];
723  M[celli] -= ((Uai[fI] & dev2GradU[fI]) & Sfb[fI])/V[celli];
724  }
725  }
726  }
727  M.correctBoundaryConditions();
728 
729  //surfaceScalarField fluxTranspose =
730  // fvc::interpolate(UaGradU, schemeName) & mesh_.Sf();
731  //forAll(mesh_.boundary(), pI)
732  //{
733  // const fvPatchScalarField& bc = nut.boundaryField()[pI];
734  // const vectorField& Sf = mesh_.boundary()[pI].Sf();
735  // if (isA<zeroGradientFvPatchScalarField>(bc))
736  // {
737  // fluxTranspose.boundaryFieldRef()[pI] =
738  // (
739  // UaGradU.boundaryField()[pI].patchInternalField()
740  // - (
741  // Ua.boundaryField()[pI].patchInternalField()
742  // & gradU_.boundaryField()[pI]
743  // )
744  // ) & Sf;
745  // }
746  // else if (isA<fixedValueFvPatchScalarField>(bc))
747  // {
748  // fluxTranspose.boundaryFieldRef()[pI] =
749  // UaGradU.boundaryField()[pI].patchInternalField() & Sf;
750  // }
751  //}
752  return
753  fvc::div(flux) - (Ua & fvc::div(snGradU)) + M;
754  //fvc::div(flux + fluxTranspose) - (Ua & fvc::div(snGradU));
755 }
756 
757 
759 (
760  const volScalarField& primalField,
761  const volScalarField& adjointField
762 ) const
763 {
764  // Grab the interpolation scheme from the primal convection term
765  tmp<surfaceInterpolationScheme<scalar>> primalInterpolationScheme
766  (
767  convectionScheme(primalField.name())
768  );
769 
770  surfaceVectorField interpolatedPrimal
771  (
772  primalInterpolationScheme().interpolate(primalField)*mesh_.Sf()
773  );
775  (
776  //reverseLinear<scalar>(mesh_).interpolate(adjointField)
777  linear<scalar>(mesh_).interpolate(adjointField)
778  *interpolatedPrimal
779  );
780 
781  const volVectorField& U = primalVars_.U();
782  forAll(mesh_.boundary(), pI)
783  {
784  const fvPatchVectorField& bc = U.boundaryField()[pI];
785  if (isA<zeroGradientFvPatchVectorField>(bc))
786  {
787  flux.boundaryFieldRef()[pI] = Zero;
788  }
789  else if (isA<fixedValueFvPatchVectorField>(bc))
790  {
791  interpolatedPrimal.boundaryFieldRef()[pI] = Zero;
792  flux.boundaryFieldRef()[pI] = Zero;
793  }
794  }
795 
796  return (-fvc::div(flux) + adjointField*fvc::div(interpolatedPrimal));
797 }
798 
799 
801 (
802  tmp<volSymmTensorField>& GbyNuMult
803 ) const
804 {
806  (
808  );
809 
810  const volVectorField& U = primalVars_.U();
811  forAll(mesh_.boundary(), pI)
812  {
813  const fvPatchVectorField& bc = U.boundaryField()[pI];
814  if (isA<zeroGradientFvPatchVectorField>(bc))
815  {
816  flux.boundaryFieldRef()[pI] = Zero;
817  }
818  else if (isA<fixedValueFvPatchVectorField>(bc))
819  {
820  flux.boundaryFieldRef()[pI] =
821  mesh_.boundary()[pI].Sf()
822  & GbyNuMult().boundaryField()[pI].patchInternalField();
823  }
824  }
825 
826  return fvc::div(flux);
827 }
828 
829 
831 (
832  tmp<volScalarField>& divUMult
833 ) const
834 {
836  (
838  );
839 
840  const volVectorField& U = primalVars_.U();
841  forAll(mesh_.boundary(), pI)
842  {
843  const fvPatchVectorField& bc = U.boundaryField()[pI];
844  if (isA<zeroGradientFvPatchVectorField>(bc))
845  {
846  flux.boundaryFieldRef()[pI] = Zero;
847  }
848  else if (isA<fixedValueFvPatchVectorField>(bc))
849  {
850  flux.boundaryFieldRef()[pI] =
851  mesh_.boundary()[pI].Sf()
852  *divUMult().boundaryField()[pI].patchInternalField();
853  }
854  }
855 
856  divUMult.clear();
857 
858  return -fvc::div(flux);
859 }
860 
861 
863 (
864  const volScalarField& primalField,
865  const volScalarField& adjointField,
866  const volScalarField& coeffField
867 ) const
868 {
869  // Note: we could grab the snGrad scheme from the diffusion term directly
870  surfaceScalarField snGradPrimal(fvc::snGrad(primalField)*mesh_.magSf());
872  (
873  reverseLinear<scalar>(mesh_).interpolate(adjointField)*snGradPrimal
874  );
875 
876  const volVectorField& U = primalVars_.U();
877  forAll(mesh_.boundary(), pI)
878  {
879  const fvPatchVectorField& bc = U.boundaryField()[pI];
880  if (!isA<coupledFvPatchVectorField>(bc))
881  {
882  flux.boundaryFieldRef()[pI] = Zero;
883  snGradPrimal.boundaryFieldRef()[pI] = Zero;
884  }
885  }
886  return (fvc::div(flux) - adjointField*fvc::div(snGradPrimal))*coeffField;
887 }
888 
889 
891 (
892  tmp<volScalarField>& mult,
893  const volScalarField& F2,
894  const volScalarField& S,
895  const volScalarField& case_1_nut,
896  const volTensorField& gradU
897 ) const
898 {
900  (
901  mult*a1_*k()*(1 - case_1_nut)/(b1_*F2*S*S*S)*twoSymm(gradU)
902  );
903  M.correctBoundaryConditions();
904  mult.clear();
905 
906  surfaceVectorField returnFieldFlux
907  (
909  );
910 
911  const volVectorField& U = primalVars_.U();
912  forAll(mesh_.boundary(), pI)
913  {
914  const fvPatchVectorField& bc = U.boundaryField()[pI];
915  if (isA<zeroGradientFvPatchVectorField>(bc))
916  {
917  returnFieldFlux.boundaryFieldRef()[pI] = Zero;
918  }
919  else if (isA<fixedValueFvPatchVectorField>(bc))
920  {
921  returnFieldFlux.boundaryFieldRef()[pI] =
922  mesh_.boundary()[pI].Sf()
923  & M.boundaryField()[pI].patchInternalField();
924  }
925  }
926 
927  // Note: potentially missing contributions form patches with a calculated
928  // nut bc
929 
930  return fvc::div(returnFieldFlux);
931 }
932 
933 
935 (
936  fvScalarMatrix& kaEqn,
937  const volScalarField& dR_dnut
938 )
939 {
940  // Add source term from the differentiation of nutWallFunction
941  scalarField& source = kaEqn.source();
943  const volScalarField& ka = this->ka();
944  const volScalarField& wa = this->wa();
945  const volScalarField& k = this->k();
946  const volScalarField& omega = this->omega();
947  forAll(nutRef().boundaryFieldRef(), patchi)
948  {
949  fvPatchScalarField& nutWall = nutRef().boundaryFieldRef()[patchi];
950  if (isA<nutkWallFunctionFvPatchScalarField>(nutWall))
951  {
952  const fvPatch& patch = mesh_.boundary()[patchi];
953  const scalarField& magSf = patch.magSf();
954 
957 
958  const scalarField& y = turbModel().y()[patchi];
959  const tmp<scalarField> tnuw = turbModel().nu(patchi);
960  const scalarField& nuw = tnuw();
961 
962  const nutWallFunctionFvPatchScalarField& nutWF =
963  refCast<nutWallFunctionFvPatchScalarField>(nutWall);
964  const wallFunctionCoefficients& wallCoeffs = nutWF.wallCoeffs();
965  const scalar Cmu = wallCoeffs.Cmu();
966  const scalar kappa = wallCoeffs.kappa();
967  const scalar E = wallCoeffs.E();
968  const scalar yPlusLam = wallCoeffs.yPlusLam();
969 
970  const scalar Cmu25 = pow025(Cmu);
971 
972  const labelList& faceCells = patch.faceCells();
973  const fvPatchScalarField& dR_dnutw =
974  dR_dnut.boundaryField()[patchi];
975  const fvPatchScalarField& omegaw = omega.boundaryField()[patchi];
976  bool addTermsFromOmegaWallFuction
977  (isA<omegaWallFunctionFvPatchScalarField>(omegaw));
978 
979  const fvPatchVectorField& Uw =
980  primalVars_.U().boundaryField()[patchi];
981  const scalarField magGradUw(mag(Uw.snGrad()));
982  forAll(nuw, facei)
983  {
984  const label celli = faceCells[facei];
985 
986  const scalar sqrtkCell(sqrt(k[celli]));
987  const scalar yPlus = Cmu25*y[facei]*sqrtkCell/nuw[facei];
988  const scalar logEyPlus = log(E*yPlus);
989  const scalar dnut_dyPlus =
990  nuw[facei]*kappa*(logEyPlus - 1)/sqr(logEyPlus);
991  const scalar dyPlus_dk =
992  Cmu25*y[facei]/(2*nuw[facei]*sqrtkCell);
993  const scalar dnut_dk = dnut_dyPlus*dyPlus_dk;
994 
995  if (yPlusLam < yPlus)
996  {
997  // Terms from nutkWallFunction
998  source[celli] -= dR_dnutw[facei]*dnut_dk*magSf[facei];
999  }
1000  if (addTermsFromOmegaWallFuction)
1001  {
1002  const scalar denom(Cmu25*kappa*y[facei]);
1003  const scalar omegaLog(sqrtkCell/denom);
1004  // Terms from omegaLog in omegaWallFunction
1005  source[celli] +=
1006  wa[celli]*omegaLog/omega[celli]
1007  /(2*sqrtkCell*denom);
1008 
1009  // Terms from G at the first cell off the wall
1010  source[celli] +=
1011  case_1_Pk_[celli]*ka[celli]*V[celli]
1012  *(
1013  (nutWall[facei] + nuw[facei])
1014  *magGradUw[facei]
1015  *Cmu25
1016  /(2.0*sqrtkCell*kappa*y[facei])
1017  );
1018 
1019  if (yPlusLam < yPlus)
1020  {
1021  source[celli] +=
1022  case_1_Pk_[celli]*ka[celli]*V[celli]
1023  *dnut_dk
1024  *magGradUw[facei]
1025  *Cmu25*sqrtkCell
1026  /(kappa*y[facei]);
1027  }
1028  }
1029  }
1030  }
1031  }
1032 }
1033 
1034 
1036 {
1038  {
1039  Info<< "Updating primal-based fields of the adjoint turbulence "
1040  << "model ..." << endl;
1041 
1042  const volVectorField& U = primalVars_.U();
1043  gradU_ = fvc::grad(U);
1044  gradOmega_ = fvc::grad(omega());
1045  gradK_ = fvc::grad(k());
1046 
1047  S2_ = 2*magSqr(symm(gradU_))
1049  S_ = sqrt(S2_);
1051 
1052  // Instead of computing G directly here, delegate to RASModelVariables
1053  // to make sure G is computed consistently when the primal fields are
1054  // averaged. The local value computed here is just used to update
1055  // the switch fields
1056  /*
1057  volScalarField G(primalVars_.turbulence()->GName(), nutRef()*GbyNu0_);
1058  omega().correctBoundaryConditions();
1059  */
1061  (
1062  IOobject
1063  (
1064  IOobject::scopedName(type(), "G"),
1065  mesh_.time().timeName(),
1066  mesh_,
1069  ),
1070  mesh_,
1072  );
1073  G.ref() = primalVars_.RASModelVariables()->G()();
1074 
1075  CDkOmega_ =
1076  (2*alphaOmega2_)*(gradK_ & gradOmega_)/omega();
1077  CDkOmegaPlus_ = max
1078  (
1079  CDkOmega_,
1080  dimensionedScalar("1.0e-10", dimless/sqr(dimTime), 1.0e-10)
1081  );
1082  F1_ = F1();
1083  F2_ = F2();
1084  case_1_Pk_ = neg(G - c1_*betaStar_*k()*omega());
1085  case_2_Pk_ = pos0(G - c1_*betaStar_*k()*omega());
1086  case_3_Pk_ = neg(a1_*omega() - b1_*F2_*S_);
1087 
1088 
1089  alphaK_ = alphaK(F1_);
1091  beta_ = beta(F1_);
1092  gamma_ = gamma(F1_);
1093 
1094  // Switch fields for F1
1095  {
1096  volScalarField arg1_C3
1097  (
1098  (scalar(1)/betaStar_)*sqrt(k())/(omega()*y_)
1099  - scalar(500)*nu()/(sqr(y_)*omega())
1100  );
1101  volScalarField arg1_C2
1102  (
1103  max
1104  (
1105  (scalar(1)/betaStar_)*sqrt(k())/(omega()*y_),
1106  scalar(500)*nu()/(sqr(y_)*omega())
1107  )
1108  - (4*alphaOmega2_)*k()/(CDkOmegaPlus_*sqr(y_))
1109  );
1110  volScalarField arg1_C1
1111  (
1112  min
1113  (
1114  max
1115  (
1116  (scalar(1)/betaStar_)*sqrt(k())/(omega()*y_),
1117  scalar(500)*nu()/(sqr(y_)*omega())
1118  ),
1119  (4*alphaOmega2_)*k()/(CDkOmegaPlus_*sqr(y_))
1120  ) - scalar(10)
1121  );
1122  volScalarField CDkOmegaPlus_C1
1123  (
1124  CDkOmega_
1125  - dimensionedScalar("1.0e-10", dimless/sqr(dimTime), 1.0e-10)
1126  );
1127 
1128  case_1_F1_ = pos(arg1_C3)*neg(arg1_C2)*neg(arg1_C1);
1129  case_2_F1_ = neg0(arg1_C3)*neg(arg1_C2)*neg(arg1_C1);
1130  case_3_F1_ = pos0(arg1_C2)*neg(arg1_C1)*pos(CDkOmegaPlus_C1);
1131  case_4_F1_ = pos0(arg1_C2)*neg(arg1_C1);
1132  }
1133 
1134  // Switch fields for nut
1135  {
1136  volScalarField nut_C1(a1_*omega() - b1_*F2_*S_);
1137  volScalarField arg2_C2
1138  (
1139  (scalar(2)/betaStar_)*sqrt(k())/(omega()*y_)
1140  - scalar(500)*nu()/(sqr(y_)*omega())
1141  );
1142  volScalarField arg2_C1
1143  (
1144  max
1145  (
1146  (scalar(2)/betaStar_)*sqrt(k())/(omega()*y_),
1147  scalar(500)*nu()/(sqr(y_)*omega())
1148  ) - scalar(100)
1149  );
1150 
1151  case_1_nut_ = pos(nut_C1);
1152  case_2_nut_ = neg0(nut_C1)*pos(arg2_C2)*neg(arg2_C1);
1153  case_3_nut_ = neg0(nut_C1)*neg0(arg2_C2)*neg(arg2_C1);
1154  }
1155 
1156  {
1157  volScalarField GPrime_C1
1158  (
1159  GbyNu0_
1160  - (c1_/a1_)*betaStar_*omega()*max(a1_*omega(), b1_*F2_*S_)
1161  );
1162  case_1_GPrime_ = neg(GPrime_C1);
1163  case_2_GPrime_ = pos0(GPrime_C1);
1164  }
1165 
1166  dnut_domega_ =
1170  DkEff_ = DkEff(F1_);
1172  changedPrimalSolution_ = false;
1173  }
1174 }
1175 
1176 
1178 (
1179  const word& varName
1180 ) const
1181 {
1183  const surfaceScalarField& phiInst = primalVars_.phiInst();
1184  word divEntry("div(" + phiInst.name() + ',' + varName +')');
1185  ITstream& divScheme = mesh_.divScheme(divEntry);
1186  // Skip the first entry which might be 'bounded' or 'Gauss'.
1187  // If it is 'bounded', skip the second entry as well
1188  word discarded(divScheme);
1189  if (discarded == "bounded")
1190  {
1191  discarded = word(divScheme);
1192  }
1194 }
1195 
1196 
1197 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1198 
1199 adjointkOmegaSST::adjointkOmegaSST
1200 (
1201  incompressibleVars& primalVars,
1202  incompressibleAdjointMeanFlowVars& adjointVars,
1203  objectiveManager& objManager,
1204  const word& adjointTurbulenceModelName,
1205  const word& modelName
1206 )
1207 :
1209  (
1210  modelName,
1211  primalVars,
1212  adjointVars,
1213  objManager,
1214  adjointTurbulenceModelName
1215  ),
1216 
1217  kappa_
1218  (
1219  dimensioned<scalar>::getOrAddToDict
1220  (
1221  "kappa",
1222  coeffDict_,
1223  0.41
1224  )
1225  ),
1226  alphaK1_
1227  (
1228  dimensioned<scalar>::getOrAddToDict
1229  (
1230  "alphaK1",
1231  this->coeffDict_,
1232  0.85
1233  )
1234  ),
1235  alphaK2_
1236  (
1237  dimensioned<scalar>::getOrAddToDict
1238  (
1239  "alphaK2",
1240  this->coeffDict_,
1241  1.0
1242  )
1243  ),
1244  alphaOmega1_
1245  (
1246  dimensioned<scalar>::getOrAddToDict
1247  (
1248  "alphaOmega1",
1249  this->coeffDict_,
1250  0.5
1251  )
1252  ),
1253  alphaOmega2_
1254  (
1255  dimensioned<scalar>::getOrAddToDict
1256  (
1257  "alphaOmega2",
1258  this->coeffDict_,
1259  0.856
1260  )
1261  ),
1262  gamma1_
1263  (
1264  dimensioned<scalar>::getOrAddToDict
1265  (
1266  "gamma1",
1267  this->coeffDict_,
1268  5.0/9.0
1269  )
1270  ),
1271  gamma2_
1272  (
1273  dimensioned<scalar>::getOrAddToDict
1274  (
1275  "gamma2",
1276  this->coeffDict_,
1277  0.44
1278  )
1279  ),
1280  beta1_
1281  (
1282  dimensioned<scalar>::getOrAddToDict
1283  (
1284  "beta1",
1285  this->coeffDict_,
1286  0.075
1287  )
1288  ),
1289  beta2_
1290  (
1291  dimensioned<scalar>::getOrAddToDict
1292  (
1293  "beta2",
1294  this->coeffDict_,
1295  0.0828
1296  )
1297  ),
1298  betaStar_
1299  (
1300  dimensioned<scalar>::getOrAddToDict
1301  (
1302  "betaStar",
1303  this->coeffDict_,
1304  0.09
1305  )
1306  ),
1307  a1_
1308  (
1309  dimensioned<scalar>::getOrAddToDict
1310  (
1311  "a1",
1312  this->coeffDict_,
1313  0.31
1314  )
1315  ),
1316  b1_
1317  (
1318  dimensioned<scalar>::getOrAddToDict
1319  (
1320  "b1",
1321  this->coeffDict_,
1322  1.0
1323  )
1324  ),
1325  c1_
1326  (
1327  dimensioned<scalar>::getOrAddToDict
1328  (
1329  "c1",
1330  this->coeffDict_,
1331  10.0
1332  )
1333  ),
1334  F3_
1335  (
1336  Switch::getOrAddToDict
1337  (
1338  "F3",
1339  this->coeffDict_,
1340  false
1341  )
1342  ),
1343 
1344  y_(primalVars_.RASModelVariables()().d()),
1345 
1346  //Primal Gradient Fields
1347  gradU_
1348  (
1349  IOobject
1350  (
1351  "rasModel::gradU",
1352  runTime_.timeName(),
1353  mesh_,
1354  IOobject::NO_READ,
1355  IOobject::NO_WRITE
1356  ),
1357  mesh_,
1359  ),
1360  gradOmega_
1361  (
1362  IOobject
1363  (
1364  "rasModel::gradOmega",
1365  runTime_.timeName(),
1366  mesh_,
1367  IOobject::NO_READ,
1368  IOobject::NO_WRITE
1369  ),
1370  mesh_,
1372  ),
1373  gradK_
1374  (
1375  IOobject
1376  (
1377  "rasModel::gradK",
1378  runTime_.timeName(),
1379  mesh_,
1380  IOobject::NO_READ,
1381  IOobject::NO_WRITE
1382  ),
1383  mesh_,
1385  ),
1386 
1387  S2_
1388  (
1389  IOobject
1390  (
1391  "S2",
1392  runTime_.timeName(),
1393  mesh_,
1394  IOobject::NO_READ,
1395  IOobject::NO_WRITE
1396  ),
1397  mesh_,
1399  ),
1400  S_
1401  (
1402  IOobject
1403  (
1404  "kOmegaSST_S",
1405  runTime_.timeName(),
1406  mesh_,
1407  IOobject::NO_READ,
1408  IOobject::NO_WRITE
1409  ),
1410  mesh_,
1412  ),
1413  GbyNu0_
1414  (
1415  IOobject
1416  (
1417  "adjointRASModel::GbyNu0",
1418  runTime_.timeName(),
1419  mesh_,
1420  IOobject::NO_READ,
1421  IOobject::NO_WRITE
1422  ),
1423  mesh_,
1425  ),
1426  CDkOmega_
1427  (
1428  IOobject
1429  (
1430  "CDkOmega_",
1431  runTime_.timeName(),
1432  mesh_,
1433  IOobject::NO_READ,
1434  IOobject::NO_WRITE
1435  ),
1436  mesh_,
1438  ),
1439  CDkOmegaPlus_
1440  (
1441  IOobject
1442  (
1443  "CDkOmegaPlus",
1444  runTime_.timeName(),
1445  mesh_,
1446  IOobject::NO_READ,
1447  IOobject::NO_WRITE
1448  ),
1449  mesh_,
1451  ),
1452  F1_
1453  (
1454  IOobject
1455  (
1456  "F1",
1457  runTime_.timeName(),
1458  mesh_,
1459  IOobject::NO_READ,
1460  IOobject::NO_WRITE
1461  ),
1462  mesh_,
1464  ),
1465  F2_
1466  (
1467  IOobject
1468  (
1469  "F2",
1470  runTime_.timeName(),
1471  mesh_,
1472  IOobject::NO_READ,
1473  IOobject::NO_WRITE
1474  ),
1475  mesh_,
1477  ),
1478  // Model Field coefficients
1479  alphaK_
1480  (
1481  IOobject
1482  (
1483  "alphaK",
1484  runTime_.timeName(),
1485  mesh_,
1486  IOobject::NO_READ,
1487  IOobject::NO_WRITE
1488  ),
1489  mesh_,
1491  ),
1492  alphaOmega_
1493  (
1494  IOobject
1495  (
1496  "alphaOmega",
1497  runTime_.timeName(),
1498  mesh_,
1499  IOobject::NO_READ,
1500  IOobject::NO_WRITE
1501  ),
1502  mesh_,
1504  ),
1505  beta_
1506  (
1507  IOobject
1508  (
1509  "beta",
1510  runTime_.timeName(),
1511  mesh_,
1512  IOobject::NO_READ,
1513  IOobject::NO_WRITE
1514  ),
1515  mesh_,
1517  ),
1518  gamma_
1519  (
1520  IOobject
1521  (
1522  "gamma",
1523  runTime_.timeName(),
1524  mesh_,
1525  IOobject::NO_READ,
1526  IOobject::NO_WRITE
1527  ),
1528  mesh_,
1530  ),
1531 
1532  case_1_F1_
1533  (
1534  IOobject
1535  (
1536  "case_1_F1",
1537  runTime_.timeName(),
1538  mesh_,
1539  IOobject::NO_READ,
1540  IOobject::NO_WRITE
1541  ),
1542  mesh_,
1544  ),
1545  case_2_F1_
1546  (
1547  IOobject
1548  (
1549  "case_2_F1",
1550  runTime_.timeName(),
1551  mesh_,
1552  IOobject::NO_READ,
1553  IOobject::NO_WRITE
1554  ),
1555  mesh_,
1557  ),
1558  case_3_F1_
1559  (
1560  IOobject
1561  (
1562  "case_3_F1",
1563  runTime_.timeName(),
1564  mesh_,
1565  IOobject::NO_READ,
1566  IOobject::NO_WRITE
1567  ),
1568  mesh_,
1570  ),
1571  case_4_F1_
1572  (
1573  IOobject
1574  (
1575  "case_4_F1",
1576  runTime_.timeName(),
1577  mesh_,
1578  IOobject::NO_READ,
1579  IOobject::NO_WRITE
1580  ),
1581  mesh_,
1583  ),
1584  case_1_Pk_
1585  (
1586  IOobject
1587  (
1588  "case_1_Pk",
1589  runTime_.timeName(),
1590  mesh_,
1591  IOobject::NO_READ,
1592  IOobject::NO_WRITE
1593  ),
1594  mesh_,
1596  ),
1597  case_2_Pk_
1598  (
1599  IOobject
1600  (
1601  "case_2_Pk",
1602  runTime_.timeName(),
1603  mesh_,
1604  IOobject::NO_READ,
1605  IOobject::NO_WRITE
1606  ),
1607  mesh_,
1609  ),
1610  case_3_Pk_
1611  (
1612  IOobject
1613  (
1614  "case_3_Pk",
1615  runTime_.timeName(),
1616  mesh_,
1617  IOobject::NO_READ,
1618  IOobject::NO_WRITE
1619  ),
1620  mesh_,
1622  ),
1623 
1624  case_1_nut_
1625  (
1626  IOobject
1627  (
1628  "case_1_nut",
1629  runTime_.timeName(),
1630  mesh_,
1631  IOobject::NO_READ,
1632  IOobject::NO_WRITE
1633  ),
1634  mesh_,
1636  ),
1637  case_2_nut_
1638  (
1639  IOobject
1640  (
1641  "case_2_nut",
1642  runTime_.timeName(),
1643  mesh_,
1644  IOobject::NO_READ,
1645  IOobject::NO_WRITE
1646  ),
1647  mesh_,
1649  ),
1650  case_3_nut_
1651  (
1652  IOobject
1653  (
1654  "case_3_nut",
1655  runTime_.timeName(),
1656  mesh_,
1657  IOobject::NO_READ,
1658  IOobject::NO_WRITE
1659  ),
1660  mesh_,
1662  ),
1663  case_1_GPrime_
1664  (
1665  IOobject
1666  (
1667  "case_1_GPrime",
1668  runTime_.timeName(),
1669  mesh_,
1670  IOobject::NO_READ,
1671  IOobject::NO_WRITE
1672  ),
1673  mesh_,
1675  ),
1676  case_2_GPrime_
1677  (
1678  IOobject
1679  (
1680  "case_2_GPrime",
1681  runTime_.timeName(),
1682  mesh_,
1683  IOobject::NO_READ,
1684  IOobject::NO_WRITE
1685  ),
1686  mesh_,
1688  ),
1689 
1690  // Zero 1rst cell field
1691  firstCellIDs_(0),
1692  zeroFirstCell_(zeroFirstCell()),
1693 
1694  // Turbulence model multipliers
1695  dnut_domega_
1696  (
1697  IOobject
1698  (
1699  "dnut_domega",
1700  runTime_.timeName(),
1701  mesh_,
1702  IOobject::NO_READ,
1703  IOobject::NO_WRITE
1704  ),
1705  mesh_,
1707  ),
1708  dnut_dk_
1709  (
1710  IOobject
1711  (
1712  "dnut_dk",
1713  runTime_.timeName(),
1714  mesh_,
1715  IOobject::NO_READ,
1716  IOobject::NO_WRITE
1717  ),
1718  mesh_,
1720  ),
1721  DOmegaEff_
1722  (
1723  IOobject
1724  (
1725  "DomegaEff",
1726  runTime_.timeName(),
1727  mesh_,
1728  IOobject::NO_READ,
1729  IOobject::NO_WRITE
1730  ),
1731  mesh_,
1732  dimensionedScalar(nutRef().dimensions(), Zero)
1733  ),
1734  DkEff_
1735  (
1736  IOobject
1737  (
1738  "DkEff",
1739  runTime_.timeName(),
1740  mesh_,
1741  IOobject::NO_READ,
1742  IOobject::NO_WRITE
1743  ),
1744  mesh_,
1745  dimensionedScalar(nutRef().dimensions(), Zero)
1746  )
1747 {
1749  adjointTMVariablesBaseNames_[0] = "ka";
1750  adjointTMVariablesBaseNames_[1] = "wa";
1751  // Read in adjoint fields
1753  (
1755  mesh_,
1756  "ka",
1757  adjointVars.solverName(),
1758  adjointVars.useSolverNameForFields()
1759  );
1761  (
1763  mesh_,
1764  "wa",
1765  adjointVars.solverName(),
1766  adjointVars.useSolverNameForFields()
1767  );
1768 
1769  setMeanFields();
1770 
1771  // No sensitivity contributions from the adjoint to the eikonal equation
1772  // for the moment
1773  includeDistance_ = false;
1774 
1775  // Update the primal related fields here so that functions computing
1776  // sensitivities have the updated fields in case of continuation
1778 }
1779 
1780 
1781 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1782 
1785  const volVectorField& Ua = adjointVars_.UaInst();
1786  return devReff(Ua);
1787 }
1788 
1789 
1791 (
1792  const volVectorField& U
1793 ) const
1794 {
1796  (
1797  IOobject
1798  (
1799  "devRhoReff",
1800  runTime_.timeName(),
1801  mesh_,
1804  ),
1806  );
1807 }
1808 
1809 
1811 {
1812  tmp<volScalarField> tnuEff = nuEff();
1813  const volScalarField& nuEff = tnuEff();
1814 
1815  return
1816  (
1817  - fvm::laplacian(nuEff, Ua)
1818  - fvc::div(nuEff*dev(fvc::grad(Ua)().T()))
1819  );
1820 
1821  /* WIP
1822  const volVectorField& U = primalVars_.U();
1823  const surfaceVectorField& Sf = mesh_.Sf();
1824  tmp<surfaceTensorField> tflux =
1825  reverseLinear<vector>(mesh_).interpolate(Ua)*Sf;
1826  surfaceTensorField& flux = tflux.ref();
1827  forAll(mesh_.boundary(), pI)
1828  {
1829  const fvPatchVectorField& Ub = U.boundaryField()[pI];
1830  if (!isA<coupledFvPatchVectorField>(Ub))
1831  {
1832  const vectorField Uai = Ua.boundaryField()[pI].patchInternalField();
1833  const vectorField& Sfb = Sf.boundaryField()[pI];
1834  flux.boundaryFieldRef()[pI] = Uai*Sfb;
1835  }
1836  }
1837  volTensorField M(nuEff*dev2(fvc::div(flux)));
1838  const DimensionedField<scalar, volMesh>& V = mesh_.V();
1839 
1840  forAll(mesh_.boundary(), pI)
1841  {
1842  const fvPatchVectorField& Ub = U.boundaryField()[pI];
1843  if (!isA<coupledFvPatchVectorField>(Ub))
1844  {
1845  const fvPatchScalarField& nuEffb = nuEff.boundaryField()[pI];
1846  const vectorField nf = mesh_.boundary()[pI].nf();
1847  const vectorField Uai = Ua.boundaryField()[pI].patchInternalField();
1848  const labelList& faceCells = mesh_.boundary()[pI].faceCells();
1849  const vectorField& Sfb = Sf.boundaryField()[pI];
1850 
1851  forAll(faceCells, fI)
1852  {
1853  const label celli = faceCells[fI];
1854  const tensor t(dev2(Uai[fI]*Sfb[fI]));
1855  M[celli] -= nuEffb[fI]*(t - nf[fI]*(nf[fI] & t))/V[celli];
1856  }
1857  }
1858  }
1859  M.correctBoundaryConditions();
1860 
1861  surfaceVectorField returnFlux =
1862  - (Sf & reverseLinear<tensor>(mesh_).interpolate(M));
1863  forAll(mesh_.boundary(), pI)
1864  {
1865  const fvPatchVectorField& Ub = U.boundaryField()[pI];
1866  if (isA<zeroGradientFvPatchVectorField>(Ub))
1867  {
1868  returnFlux.boundaryFieldRef()[pI] = Zero;
1869  }
1870  else if (isA<fixedValueFvPatchVectorField>(Ub))
1871  {
1872  const scalarField& deltaCoeffs = mesh_.boundary()[pI].deltaCoeffs();
1873  const fvPatchScalarField& nuEffb = nuEff.boundaryField()[pI];
1874  const vectorField Uai = Ua.boundaryField()[pI].patchInternalField();
1875  const vectorField nf = mesh_.boundary()[pI].nf();
1876  const vectorField& Sfb = Sf.boundaryField()[pI];
1877 
1878  returnFlux.boundaryFieldRef()[pI] =
1879  - (Sfb & M.boundaryField()[pI].patchInternalField())
1880  + nuEffb*deltaCoeffs*(nf & dev2(Uai*Sfb));
1881  }
1882  }
1883 
1884  return
1885  (
1886  - fvm::laplacian(nuEff, Ua)
1887  + fvc::div(returnFlux)
1888  );
1889  */
1890 }
1891 
1894 {
1895  return (ka()*gradK_ + wa()*gradOmega_);
1896 }
1897 
1898 
1900 {
1901  tmp<volVectorField> tmeanFlowSource
1902  (
1904  (
1905  IOobject
1906  (
1907  "adjointMeanFlowSource" + type(),
1908  mesh_.time().timeName(),
1909  mesh_,
1912  ),
1913  mesh_,
1915  )
1916  );
1917  volVectorField& meanFlowSource = tmeanFlowSource.ref();
1918 
1919  // Contributions from the convection terms of the turbulence model
1920  meanFlowSource +=
1922  + convectionMeanFlowSource(k(), ka());
1923 
1924  // Contributions from GbyNu, including gradU
1925  tmp<volSymmTensorField> twoSymmGradU(twoSymm(gradU_));
1926  tmp<volSymmTensorField> GbyNuMult
1927  (
1928  // First part of GPrime and G from Pk
1929  2.*dev(twoSymmGradU())*zeroFirstCell_
1931  // Second part of GPrime
1932  + twoSymmGradU()*wa()*zeroFirstCell_*gamma_*case_2_GPrime_
1934  );
1935  twoSymmGradU.clear();
1936  meanFlowSource += GMeanFlowSource(GbyNuMult);
1937  GbyNuMult.clear();
1938 
1939  // Contributions from divU
1940  tmp<volScalarField> divUMult
1941  (
1942  (2.0/3.0)*(zeroFirstCell_*wa()*omega()*gamma_ + ka()*k())
1943  );
1944  meanFlowSource += divUMeanFlowSource(divUMult);
1945 
1946  // Contributions from S2, existing in nut
1947  const volVectorField& U = primalVars_.U();
1948  const volVectorField& Ua = adjointVars_.UaInst();
1949  tmp<volScalarField> nutMeanFlowSourceMult
1950  (
1951  // nut in the diffusion coefficients
1954  + dNutdbMult(U, Ua, nutRef(), "coeffsDiff")
1955  // nut in G
1957  );
1958  meanFlowSource +=
1960  (
1961  nutMeanFlowSourceMult,
1962  F2_,
1963  S_,
1964  case_1_nut_,
1965  gradU_
1966  );
1967 
1968  // G at the first cell includes mag(U.snGrad())
1969  // Add term here
1970  forAll(omega().boundaryFieldRef(), patchi)
1971  {
1972  fvPatchScalarField& omegaWall = omega().boundaryFieldRef()[patchi];
1973  if (isA<omegaWallFunctionFvPatchScalarField>(omegaWall))
1974  {
1975  const fvPatch& patch = mesh_.boundary()[patchi];
1976 
1977  const autoPtr<incompressible::turbulenceModel>& turbModel =
1979 
1980  const scalarField& y = turbModel().y()[patchi];
1981  const tmp<scalarField> tnuw = turbModel().nu(patchi);
1982  const scalarField& nuw = tnuw();
1983 
1984  const nutWallFunctionFvPatchScalarField& nutw =
1985  refCast<nutWallFunctionFvPatchScalarField>
1986  (nutRef().boundaryFieldRef()[patchi]);
1987  const wallFunctionCoefficients& wallCoeffs = nutw.wallCoeffs();
1988  const scalar Cmu = wallCoeffs.Cmu();
1989  const scalar kappa = wallCoeffs.kappa();
1990  const scalar Cmu25 = pow025(Cmu);
1991 
1992  const labelList& faceCells = patch.faceCells();
1993 
1994  const fvPatchVectorField& Uw =
1995  primalVars_.U().boundaryField()[patchi];
1996  vectorField snGradUw(Uw.snGrad());
1997  const scalarField& deltaCoeffs = patch.deltaCoeffs();
1998  forAll(faceCells, facei)
1999  {
2000  const label celli = faceCells[facei];
2001  // Volume will be added when meanFlowSource is added to UaEqn
2002  meanFlowSource[celli] +=
2003  ka()[celli]*case_1_Pk_[celli]
2004  *(nutw[facei] + nuw[facei])
2005  *snGradUw[facei].normalise()
2006  *Cmu25*sqrt(k()[celli])
2007  *deltaCoeffs[facei]
2008  /(kappa*y[facei]);
2009  }
2010  }
2011  }
2012  return tmeanFlowSource;
2013 }
2014 
2015 
2016 tmp<volScalarField> adjointkOmegaSST::nutJacobianTMVar1() const
2017 {
2018  // Compute dnut_dk anew since the current copy of dnut_dk_
2019  // might not be updated
2020  const volVectorField& U = primalVars_.U();
2021  tmp<volScalarField> S2
2022  (
2023  2*magSqr(symm(fvc::grad(U)))
2025  );
2026  volScalarField S(sqrt(S2));
2027  volScalarField F2(this->F2());
2028 
2029  // Computation of nut switches
2030  volScalarField nut_C1(a1_*omega() - b1_*F2*S);
2031  volScalarField arg2_C2
2032  (
2033  (scalar(2)/betaStar_)*sqrt(k())/(omega()*y_)
2034  - scalar(500)*nu()/(sqr(y_)*omega())
2035  );
2036  volScalarField arg2_C1
2037  (
2038  max
2039  (
2040  (scalar(2)/betaStar_)*sqrt(k())/(omega()*y_),
2041  scalar(500)*nu()/(sqr(y_)*omega())
2042  ) - scalar(100)
2043  );
2044  volScalarField case_2_nut(neg0(nut_C1)*pos(arg2_C2)*neg(arg2_C1));
2045 
2046  return dnut_dk(F2, S, case_2_nut);
2047 }
2048 
2049 
2051 {
2052  // Compute dnut_omega anew since the current copy of dnut_domega_
2053  // might not be updated
2054  const volVectorField& U = primalVars_.U();
2056  (
2057  2*magSqr(symm(fvc::grad(U)))
2059  );
2060  volScalarField S(sqrt(S2));
2061  volScalarField F2(this->F2());
2062 
2063  // Computation of nut switches
2064  volScalarField nut_C1(a1_*omega() - b1_*F2*S);
2065  volScalarField arg2_C2
2066  (
2067  (scalar(2)/betaStar_)*sqrt(k())/(omega()*y_)
2068  - scalar(500)*nu()/(sqr(y_)*omega())
2069  );
2070  volScalarField arg2_C1
2071  (
2072  max
2073  (
2074  (scalar(2)/betaStar_)*sqrt(k())/(omega()*y_),
2075  scalar(500)*nu()/(sqr(y_)*omega())
2076  ) - scalar(100)
2077  );
2078  volScalarField case_1_nut(pos(nut_C1));
2079  volScalarField case_2_nut(neg0(nut_C1)*pos(arg2_C2)*neg(arg2_C1));
2080  volScalarField case_3_nut(neg0(nut_C1)*neg0(arg2_C2)*neg(arg2_C1));
2081 
2082  return dnut_domega(F2, S, case_1_nut, case_2_nut, case_3_nut);
2083 }
2084 
2085 
2087 (
2088  tmp<volScalarField>& dNutdUMult
2089 ) const
2090 {
2091  const volVectorField& U = primalVars_.U();
2092  volTensorField gradU(fvc::grad(U));
2094  (
2095  2*magSqr(symm(gradU))
2097  );
2098  volScalarField S(sqrt(S2));
2099  volScalarField F2(this->F2());
2100 
2101  // Computation of nut switches
2102  volScalarField nut_C1(a1_*omega() - b1_*F2*S);
2103  volScalarField case_1_nut(pos(nut_C1));
2104 
2105  return nutMeanFlowSource(dNutdUMult, F2, S, case_1_nut, gradU);
2106 }
2107 
2108 
2110 {
2111  return
2112  (
2113  alphaK_.boundaryField()[patchI]*nutRef().boundaryField()[patchI]
2114  + nu()().boundaryField()[patchI]
2115  );
2116 }
2117 
2118 
2120 {
2121  return
2122  (
2123  alphaOmega_.boundaryField()[patchI]*nutRef().boundaryField()[patchI]
2124  + nu()().boundaryField()[patchI]
2125  );
2126 }
2127 
2128 
2130 {
2132 
2133  if (!adjointTurbulence_)
2134  {
2135  return;
2136  }
2137 
2139 
2140  // Primal and adjoint fields
2141  const volVectorField& U = primalVars_.U();
2143  volScalarField dR_dnut(this->dR_dnut());
2145 
2147 
2148  tmp<fvScalarMatrix> waEqn
2149  (
2150  fvm::div(-phi, wa())
2154  + waEqnSourceFromF1()
2156  + fvm::Sp(zeroFirstCell_*scalar(2.)*beta_*omega(), wa())
2157  + fvm::SuSp
2158  (
2159  zeroFirstCell_()*gamma_*((2.0/3.0)*divU - dGPrime_domega().ref()()),
2160  wa()
2161  )
2162  - (case_2_Pk_*c1_ - scalar(1))*betaStar_*k()*ka()
2163  ==
2164  fvOptions(wa())
2165  );
2166 
2167  // Boundary manipulate changes the diagonal component, so relax has to
2168  // come after that
2169  waEqn.ref().boundaryManipulate(wa().boundaryFieldRef());
2170  waEqn.ref().relax();
2171 
2172  // Sources from the objective should be added after the boundary
2173  // manipulation
2174  objectiveManager_.addSource(waEqn.ref());
2175  fvOptions.constrain(waEqn.ref());
2176  waEqn.ref().solve();
2177 
2178  // Adjoint Turbulent kinetic energy equation
2179  tmp<fvScalarMatrix> kaEqn
2180  (
2181  fvm::div(-phi, ka())
2182  + fvm::SuSp(fvc::div(phi), ka())
2183  - fvm::laplacian(DkEff_, ka())
2184  + fvm::Sp(betaStar_*omega(), ka())
2185  - case_2_Pk_()*c1_*betaStar_*omega()()*ka()()
2186  + fvm::SuSp(scalar(2.0/3.0)*divU, ka())
2188  + kaEqnSourceFromF1()
2189  + dR_dnut()*dnut_dk_()
2190  - zeroFirstCell_()*gamma_*dGPrime_dk().ref()()*wa()
2191  ==
2192  fvOptions(ka())
2193  );
2194 
2195  kaEqn.ref().relax();
2196  kaEqn.ref().boundaryManipulate(ka().boundaryFieldRef());
2197  addWallFunctionTerms(kaEqn.ref(), dR_dnut);
2198  fvOptions.constrain(kaEqn.ref());
2199  // Add sources from the objective functions
2200  objectiveManager_.addSource(kaEqn.ref());
2201 
2202  kaEqn.ref().solve();
2203 
2205  {
2206  dimensionedScalar maxwa = max(mag(wa()));
2208  Info<< "Max mag (" << wa().name() << ") = " << maxwa.value() << endl;
2209  Info<< "Max mag (" << ka().name() << ") = " << maxka.value() << endl;
2210  }
2211 }
2212 
2215 {
2216  return adjMomentumBCSourcePtr_();
2217 }
2218 
2219 
2221 {
2224 
2225  forAll(mesh_.boundary(), patchi)
2226  {
2227  vectorField nf(mesh_.boundary()[patchi].nf());
2228  wallShapeSens[patchi] = nf & FITerm.boundaryField()[patchi];
2229  }
2230  return wallShapeSens;
2231 }
2232 
2235 {
2236  return wallFloCoSensitivitiesPtr_();
2237 }
2238 
2239 
2241 {
2243  (
2244  IOobject
2245  (
2246  "adjointEikonalSource" + type(),
2247  runTime_.timeName(),
2248  mesh_,
2251  ),
2252  mesh_,
2254  );
2255 }
2256 
2257 
2259 {
2260  const volVectorField& U = primalVars_.U();
2261  const volScalarField& kInst =
2262  primalVars_.RASModelVariables()->TMVar1Inst();
2263  const volScalarField& omegaInst =
2264  primalVars_.RASModelVariables()->TMVar2Inst();
2265 
2266  tmp<volScalarField> arg1 = min
2267  (
2268  min
2269  (
2270  max
2271  (
2272  (scalar(1)/betaStar_)*sqrt(k())/(omega()*y_),
2273  scalar(500)*nu()/(sqr(y_)*omega())
2274  ),
2275  (4*alphaOmega2_)*k()/(CDkOmegaPlus_*sqr(y_))
2276  ),
2277  scalar(10)
2278  );
2279 
2280  // Interpolation schemes used by the primal convection terms
2281  auto kScheme(convectionScheme(kInst.name()));
2282  auto omegaScheme(convectionScheme(omegaInst.name()));
2283  const surfaceVectorField& Sf = mesh_.Sf();
2284  tmp<volTensorField> tFISens
2285  (
2287  (
2288  IOobject
2289  (
2290  type() + "FISensTerm",
2291  mesh_.time().timeName(),
2292  mesh_,
2295  ),
2296  mesh_,
2299  )
2300  );
2301  volTensorField& FISens = tFISens.ref();
2302  FISens =
2303  // k convection
2304  - ka()*fvc::div
2305  (
2306  kScheme().interpolate(k())
2308  )
2309  // k diffusion
2310  + ka()*T(fvc::grad(DkEff_*gradK_))
2311  - DkEff_*(fvc::grad(ka())*gradK_)
2312  // omega convection
2314  (
2315  omegaScheme().interpolate(omega())
2317  )
2318  // omega diffusion
2321  // terms including GbyNu0
2322  + (
2324  + case_1_Pk_*ka()*nutRef()
2326  // S2 (includes contribution from nut in UEqn as well)
2327  + (
2328  dR_dnut()*a1_*k()/(b1_*S_*S_*S_*F2_)
2330  *(c1_/a1_)*betaStar_*omega()*b1_*F2_/S_
2331  )*T(gradU_ & twoSymm(gradU_))*(1. - case_1_nut_)
2332  // CDkOmega in omegaEqn
2333  + 2.*wa()*(1. - F1_)*alphaOmega2_/omega()*zeroFirstCell_
2335  // F1
2336  - dR_dF1()
2337  *(dF1_dGradK(arg1)*gradK_ + dF1_dGradOmega(arg1)*gradOmega_);
2338 
2340 
2341  return tFISens;
2342 }
2343 
2344 
2346 (
2347  const word& designVarsName
2348 ) const
2349 {
2351  auto tres(tmp<scalarField>::New(mesh_.nCells(), Zero));
2352 
2353  // Sensitivity from the source term in the k equation
2354  scalarField auxSens
2355  (k().primitiveField()*ka().primitiveField());
2357  (
2358  tres.ref(), auxSens, fvOptions, k().name(), designVarsName
2359  );
2360 
2361  // Sensitivity from the source term in the omega equation
2362  auxSens = omega().primitiveField()*wa().primitiveField();
2364  (
2365  tres.ref(), auxSens, fvOptions, omega().name(), designVarsName
2366  );
2367 
2368  return tres;
2369 }
2370 
2371 
2373 {
2376 }
2377 
2378 
2380 {
2381  if (adjointRASModel::read())
2382  {
2384  alphaK1_.readIfPresent(this->coeffDict());
2385  alphaK2_.readIfPresent(this->coeffDict());
2388  gamma1_.readIfPresent(this->coeffDict());
2389  gamma2_.readIfPresent(this->coeffDict());
2390  beta1_.readIfPresent(this->coeffDict());
2391  beta2_.readIfPresent(this->coeffDict());
2393  a1_.readIfPresent(this->coeffDict());
2394  b1_.readIfPresent(this->coeffDict());
2395  c1_.readIfPresent(this->coeffDict());
2396  F3_.readIfPresent("F3", this->coeffDict());
2397 
2398  return true;
2399  }
2400  else
2401  {
2402  return false;
2403  }
2404 }
2405 
2406 
2407 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
2408 
2409 } // End namespace adjointRASModels
2410 } // End namespace incompressible
2411 } // End namespace Foam
2412 
2413 // ************************************************************************* //
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:153
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:76
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:72
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:320
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.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
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:1095
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...
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:1261
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:180
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)
Do not request registration (bool: false)
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