Loading...
Searching...
No Matches
adjointSpalartAllmaras.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) 2007-2023 PCOpt/NTUA
9 Copyright (C) 2013-2023 FOSS GP
10 Copyright (C) 2019-2023 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
32#include "wallDist.H"
33#include "wallFvPatch.H"
36#include "coupledFvPatch.H"
37#include "ATCModel.H"
38#include "fvOptions.H"
39#include "sensitivityTopO.H"
40
41// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42
43namespace Foam
44{
47namespace adjointRASModels
48{
49
50// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
51
54(
58);
59
60
61// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
62
63// * * * * * * * * * * * * Primal Spalart - Allmaras * * * * * * * * * * * * //
69
70
79(
80 const volScalarField& chi,
82) const
83{
84 return 1.0 - chi/(1.0 + chi*fv1);
85}
86
87
89(
90 const volScalarField& chi,
91 const volScalarField& fv1
92) const
93{
94 volScalarField Omega(::sqrt(2.0)*mag(skew(gradU_)));
95
96 return
97 (
98 max
99 (
100 Omega
102 Cs_*Omega
103 )
104 );
105}
106
107
109(
110 const volScalarField& Stilda
111) const
112{
113 auto tr =
115 (
116 min
117 (
119 scalar(10)
120 )
121 );
122 tr.ref().boundaryFieldRef() == Zero;
123
124 return tr;
125}
126
127
129(
130 const volScalarField& Stilda
131) const
133 const volScalarField g(r_ + Cw2_*(pow6(r_) - r_));
134
135 return g*pow((1.0 + pow6(Cw3_))/(pow6(g) + pow6(Cw3_)), 1.0/6.0);
136}
137
138
141 return
143 ("DnuTildaEff", (nuTilda() + this->nu())/sigmaNut_);
144}
145
148{
149 return primalVars_.RASModelVariables()().TMVar1();
150}
151
152
155 return primalVars_.RASModelVariables()().nutRef();
156}
157
158
159// * * * * * * * * * * * Adjoint Spalart - Allmaras * * * * * * * * * * * * //
160
162(
163 const volScalarField& chi
164) const
166 volScalarField chi3(pow3(chi));
167
168 return 3.0*pow3(Cv1_)*sqr(chi/(chi3 + pow3(Cv1_)));
169}
170
171
173(
174 const volScalarField& chi,
175 const volScalarField& fv1,
176 const volScalarField& dFv1dChi
177) const
178{
179 return (chi*chi*dFv1dChi - 1.)/sqr(1. + chi*fv1);
180}
181
182
184(
185 const volScalarField& Omega,
186 const volScalarField& fv2
187) const
188{
189 volScalarField fieldSwitch
190 (
191 Omega + fv2*nuTilda()/sqr(kappa_*y_) - Cs_*Omega
192 );
193
194 return pos(fieldSwitch) + neg(fieldSwitch)*Cs_;
195}
196
197
199(
200 const volScalarField& Omega,
201 const volScalarField& fv2,
202 const volScalarField& dFv2dChi
203) const
204{
206 volScalarField fieldSwitch(Omega + fv2*nuTilda()*invDenom - Cs_*Omega);
207
208 return pos(fieldSwitch)*(dFv2dChi*nuTilda()*invDenom/nu() + fv2*invDenom);
209}
210
211
213(
214 const volScalarField& Omega,
215 const volScalarField& fv2
216) const
217{
219 volScalarField fieldSwitch(Omega + aux - Cs_*Omega);
220
221 return - 2.*pos(fieldSwitch)*aux/y_;
222}
223
224
226(
227 const volScalarField& Stilda
228) const
229{
230 tmp<volScalarField> tdrdNutilda
231 (
233 *(scalar(10) - r_)/(scalar(10) - r_ + SMALL)
234 );
235 tdrdNutilda.ref().boundaryFieldRef() == Zero;
236
237 return tdrdNutilda;
238}
239
240
242(
243 const volScalarField& Stilda
244) const
245{
246 tmp<volScalarField> tdrdStilda
247 (
249 *(scalar(10) - r_)/(scalar(10) - r_ + SMALL)
250 );
251 tdrdStilda.ref().boundaryFieldRef() == Zero;
252
253 return tdrdStilda;
254}
255
256
258(
259 const volScalarField& Stilda
260) const
261{
262 tmp<volScalarField> tdrdDelta
263 (
265 *(scalar(10) - r_)/(scalar(10) - r_ + SMALL)
266 );
267 tdrdDelta.ref().boundaryFieldRef() == Zero;
268
269 return tdrdDelta;
270}
271
272
274(
275 const volScalarField& Stilda
276) const
277{
278 volScalarField g(r_ + Cw2_*(pow6(r_) - r_));
279
280 dimensionedScalar pow6Cw3 = pow6(Cw3_);
281 volScalarField pow6g(pow6(g));
282
283 return
284 pow6Cw3/(pow6g + pow6Cw3)
285 *pow((1.0 + pow6Cw3)/(pow6g + pow6Cw3), 1.0/6.0)
286 *(1.0 + Cw2_*(6.0*pow5(r_) - 1.0));
287}
288
289
291(
292 const volScalarField& Stilda,
293 const volScalarField& dfwdr,
294 const volScalarField& dStildadNuTilda
295) const
296{
298
299 return
300 dfwdr*(dr_dNuTilda(Stilda) + dr_dStilda(Stilda)*dStildadNuTilda);
301}
302
303
305(
306 const volScalarField& Stilda,
307 const volScalarField& dfwdr,
308 const volScalarField& dStildadOmega
309) const
310{
311 return dfwdr*dr_dStilda(Stilda)*dStildadOmega;
312}
313
314
316(
317 const volScalarField& Stilda,
318 const volScalarField& dfwdr,
319 const volScalarField& dStildadDelta
320) const
321{
322 return dfwdr*(dr_dDelta(Stilda) + dr_dStilda(Stilda)*dStildadDelta);
323}
324
325
327(
328 const volScalarField& fw,
329 const volScalarField& dfwdNuTilda
330) const
331{
332 return Cw1_*(nuTilda()*dfwdNuTilda + fw)/sqr(y_);
333}
334
335
337(
338 const volScalarField& dStildadNuTilda
339) const
340{
341 return - Cb1_*dStildadNuTilda;
342}
343
344
346(
347 const volScalarField& fv1,
348 const volScalarField& dFv1dChi
349) const
350{
351 return dFv1dChi*nuTilda()/nu() + fv1;
352}
353
354
356{
357 // Store boundary field of the conservative part,
358 // for use in adjoint outlet boundary conditions
359 forAll(mesh_.boundary(), pI)
360 {
361 const fvPatch& patch = mesh_.boundary()[pI];
362 if(!isA<coupledFvPatch>(patch))
363 {
364 tmp<vectorField> tnf = patch.nf();
367 *nuaTilda().boundaryField()[pI];
369 }
370
372}
373
374
376{
378 {
379 Info<< "Updating primal-based fields of the adjoint turbulence "
380 << "model ..." << endl;
381
382 // Grab references
383 const volVectorField& U = primalVars_.U();
384
385 // Update gradient fields
386 gradU_ = mask_*fvc::grad(U, "gradUStilda");
388
389 const volScalarField Omega(::sqrt(2.0)*mag(skew(gradU_)));
390
391 // Primal SA fields
392 volScalarField chi(this->chi());
393 volScalarField fv1(this->fv1(chi));
394 volScalarField fv2(this->fv2(chi, fv1));
395 Stilda_ = Stilda(chi, fv1);
396 r_ = r(Stilda_);
397 fw_ = this->fw(Stilda_);
398
399 // Derivatives of primal fields wrt to nuTilda
403 (this->dStilda_dNuTilda(Omega, fv2, dFv2_dChi));
407
408 // Fields to be used in the nuaTilda equation
410 symm(mask_*fvc::grad(U, "adjointProductionU"));
411
413 nuTilda()
414 *(
417 );
418
420
421 // Constant multiplier in the adjoint momentum source term
425
427 2.*skew(gradU_)
428 /(Omega + dimensionedScalar("SMALL", Omega.dimensions(), SMALL))
429 *(
432 );
433
434 // Set changedPrimalSolution_ to false to avoid recomputing these
435 // fields unless the primal has changed
437 }
438}
439
440
442{
444 {
446 }
447
449 (
450 "unitMask",
452 mesh_,
453 dimensionedScalar("unit", dimless, scalar(1))
454 );
455}
456
457
458// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
459
460adjointSpalartAllmaras::adjointSpalartAllmaras
461(
462 incompressibleVars& primalVars,
464 objectiveManager& objManager,
465 const word& adjointTurbulenceModelName,
466 const word& modelName
467)
468:
470 (
471 modelName,
472 primalVars,
473 adjointVars,
474 objManager,
475 adjointTurbulenceModelName
476 ),
477
478 sigmaNut_
479 (
480 dimensioned<scalar>::getOrAddToDict
481 (
482 "sigmaNut",
483 this->coeffDict_,
484 0.66666
485 )
486 ),
487 kappa_
488 (
489 dimensioned<scalar>::getOrAddToDict
490 (
491 "kappa",
492 this->coeffDict_,
493 0.41
494 )
495 ),
496
497 Cb1_
498 (
499 dimensioned<scalar>::getOrAddToDict
500 (
501 "Cb1",
502 this->coeffDict_,
503 0.1355
504 )
505 ),
506 Cb2_
507 (
508 dimensioned<scalar>::getOrAddToDict
509 (
510 "Cb2",
511 this->coeffDict_,
512 0.622
513 )
514 ),
515 Cw1_(Cb1_/sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_),
516 Cw2_
517 (
518 dimensioned<scalar>::getOrAddToDict
519 (
520 "Cw2",
521 this->coeffDict_,
522 0.3
523 )
524 ),
525 Cw3_
526 (
527 dimensioned<scalar>::getOrAddToDict
528 (
529 "Cw3",
530 this->coeffDict_,
531 2.0
532 )
533 ),
534 Cv1_
535 (
536 dimensioned<scalar>::getOrAddToDict
537 (
538 "Cv1",
539 this->coeffDict_,
540 7.1
541 )
542 ),
543 Cs_
544 (
545 dimensioned<scalar>::getOrAddToDict
546 (
547 "Cs",
548 this->coeffDict_,
549 0.3
550 )
551 ),
552
553 limitAdjointProduction_
554 (
555 coeffDict_.getOrDefault("limitAdjointProduction", true)
556 ),
557
558 y_(primalVars_.RASModelVariables()().d()),
559
560 mask_(allocateMask()),
561
562 symmAdjointProductionU_
563 (
565 (
566 "symmAdjointProductionU",
567 runTime_.timeName(),
568 mesh_,
569 IOobject::NO_READ,
570 IOobject::NO_WRITE
571 ),
572 mesh_,
574 ),
575
576 productionDestructionSource_
577 (
579 (
580 "productionDestructionSource",
581 runTime_.timeName(),
582 mesh_,
583 IOobject::NO_READ,
584 IOobject::NO_WRITE
585 ),
586 mesh_,
588 ),
589
590 Stilda_
591 (
593 (
594 "Stilda",
595 runTime_.timeName(),
596 mesh_,
597 IOobject::NO_READ,
598 IOobject::NO_WRITE
599 ),
600 mesh_,
602 ),
603
604 r_
605 (
607 (
608 "r",
609 runTime_.timeName(),
610 mesh_,
611 IOobject::NO_READ,
612 IOobject::NO_WRITE
613 ),
614 mesh_,
616 ),
617
618 fw_
619 (
621 (
622 "fw",
623 runTime_.timeName(),
624 mesh_,
625 IOobject::NO_READ,
626 IOobject::NO_WRITE
627 ),
628 mesh_,
630 ),
631
632 Cdnut_
633 (
635 (
636 "Cdnut",
637 runTime_.timeName(),
638 mesh_,
639 IOobject::NO_READ,
640 IOobject::NO_WRITE
641 ),
642 mesh_,
644 ),
645
646 momentumSourceMult_
647 (
649 (
650 "momentumSourceMult",
651 runTime_.timeName(),
652 mesh_,
653 IOobject::NO_READ,
654 IOobject::NO_WRITE
655 ),
656 mesh_,
658 ),
659
660 gradU_(fvc::grad(primalVars.U())),
661 gradNuTilda_(fvc::grad(nuTilda())),
662 minStilda_("SMALL", Stilda_.dimensions(), SMALL)
663{
665 adjointTMVariablesBaseNames_[0] = "nuaTilda";
666 // Read nuaTilda field and reset pointer to the first
667 // adjoint turbulence model variable
669 (
671 mesh_,
672 "nuaTilda",
673 adjointVars.solverName(),
674 adjointVars.useSolverNameForFields()
675 );
676
678
679 // Set the includeDistance to true, to allow for the automatic solution
680 // of the adjoint eikonal equation when computing sensitivities
681 includeDistance_ = true;
682
683 // Update the primal related fields here so that functions computing
684 // sensitivities have the updated fields in case of continuation
686}
687
688
689// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
690
692{
693 const volVectorField& Ua = adjointVars_.UaInst();
694 return devReff(Ua);
695}
696
697
699(
700 const volVectorField& U
701) const
702{
704 (
705 "devRhoReff",
708 );
709}
710
711
713{
714 tmp<volScalarField> tnuEff(nuEff());
715 const volScalarField& nuEff = tnuEff();
716
717 return
719 - fvm::laplacian(nuEff, Ua)
720 - fvc::div(nuEff*dev(fvc::grad(Ua)().T()))
721 );
722}
723
724
726{
728 (
730 "adjointSpalartAllmaras::addMomentumSource"
731 );
732 // cm formulation
733 //return (- nuTilda()*fvc::grad(nuaTilda() - conservativeMomentumSource());
734
735 // ncm formulation
737}
738
739
741{
742 volScalarField chi(this->chi());
743 volScalarField fv1(this->fv1(chi));
745
746 return dnut_dNuTilda(fv1, dFv1_dChi);
747}
748
749
751{
752 auto tdiffCoeff =
754
755 scalarField& diffCoeff = tdiffCoeff.ref();
756
757 diffCoeff =
758 (nuTilda().boundaryField()[patchI] + nu()().boundaryField()[patchI])
759 /sigmaNut_.value();
760
761 return tdiffCoeff;
762}
763
764
767{
768 // Computed in conservativeMomentumSource
770}
771
772
774{
776
777 forAll(mesh_.boundary(), patchI)
778 {
779 const fvPatch& patch = mesh_.boundary()[patchI];
780
781 tmp<vectorField> tnf = patch.nf();
782 if (isA<wallFvPatch>(patch) && patch.size() != 0)
783 {
784 wallShapeSens[patchI] =
785 - nuaTilda().boundaryField()[patchI].snGrad()
786 *diffusionCoeffVar1(patchI)
787 *nuTilda().boundaryField()[patchI].snGrad()*tnf;
789 }
790
791 return wallShapeSens;
792}
793
794
796{
798
799 forAll(mesh_.boundary(), patchI)
800 {
801 tmp<vectorField> tnf = mesh_.boundary()[patchI].nf();
802
803 wallFloCoSens[patchI] =
804 nuaTilda().boundaryField()[patchI]
805 *nuTilda().boundaryField()[patchI]*tnf;
806 }
807
808 return wallFloCoSens;
809}
810
811
813{
814 const volVectorField& U = primalVars_.U();
815 const volVectorField& Ua = adjointVars_.Ua();
816
817 // Primal SA fields
818 volScalarField chi(this->chi());
819 volScalarField fv1(this->fv1(chi));
820 volScalarField fv2(this->fv2(chi, fv1));
821 volScalarField Omega(::sqrt(2.0)*mag(gradU_));
822
823 // Derivatives of primal fields wrt to nuTilda
829 (this->dfw_dDelta(Stilda_, dfw_dr, dStilda_dDelta));
830
831 auto tadjointEikonalSource =
833 (
834 "adjointEikonalSource" + type(),
835 (
837 + Cw1_*sqr(nuTilda()/y_)*(dfw_dDelta - 2.*fw_/y_)
838 )*nuaTilda()
839 );
840 volScalarField& adjointEikonalSource = tadjointEikonalSource.ref();
841
842 // if wall functions are used, add appropriate source terms
844 SAwallFunctionPatchField;
845
846 const volScalarField::Boundary& nutBoundary = nutRef().boundaryField();
847 const scalarField& V = mesh_.V().field();
848
849 tmp<volScalarField> tnuEff = nuEff();
850 const volScalarField& nuEff = tnuEff();
851
852 forAll(nutBoundary, patchi)
853 {
854 const fvPatch& patch = mesh_.boundary()[patchi];
855 if
856 (
857 isA<SAwallFunctionPatchField>(nutBoundary[patchi])
858 && patch.size() != 0
859 )
860 {
861 const scalar kappa_(0.41);
862 const scalar E_(9.8);
863 tmp<vectorField> tnf = patch.nf();
864 const vectorField& nf = tnf.ref();
865 const scalarField& magSf = patch.magSf();
866
867 const fvPatchVectorField& Up = U.boundaryField()[patchi];
868 const fvPatchVectorField& Uap = Ua.boundaryField()[patchi];
869 const vectorField Uc(Up.patchInternalField());
870 const vectorField Uc_t(Uc - (Uc & nf)*nf);
871
872 // By convention, tf has the direction of the tangent
873 // PRIMAL velocity at the first cell off the wall
874 const vectorField tf(Uc_t/mag(Uc_t));
875
876 const scalarField nuw(nuEff.boundaryField()[patchi]);
877 const scalarField nu(this->nu()().boundaryField()[patchi]);
878 const fvPatchScalarField& yC = y()[patchi];
879
880 const scalarField magGradU(mag(Up.snGrad()));
881
882 // Note: What happens in separation?? sign change needed
883 const scalarField vtau(sqrt(nuw*magGradU));
884
885 // Note: mag for positive uPlus
886 const scalarField uPlus(mag(Uc)/vtau);
887
888 const scalarField yPlus(yC*vtau/nu);
889 const scalarField kUu(min(kappa_*uPlus, scalar(50)));
890 const scalarField auxA
891 ((kappa_/E_)*(exp(kUu) - 1 - kUu - 0.5*kUu*kUu));
892 const scalarField Cwf_d(sqr(vtau)/nu/(yPlus+uPlus*(1 + auxA)));
893
894 // Tangential components are according to tf
896 (
898 (
899 "objectiveManager" + objectiveManager_.adjointSolverName(),
901 "incompressible",
902 patch
903 )
904 );
905 tmp<vectorField> tsource(boundaryContrPtr->normalVelocitySource());
906
907 const scalarField rt(tsource() & tf);
908 const scalarField Uap_t(Uap & tf);
909
910 const labelUList& faceCells = patch.faceCells();
911 forAll(faceCells, faceI)
912 {
913 label cellI = faceCells[faceI];
914 adjointEikonalSource[cellI] -=
915 2.*( rt[faceI] + Uap_t[faceI] )
916 *vtau[faceI]*Cwf_d[faceI]*magSf[faceI]
917 /V[cellI]; // Divide with cell volume since the term
918 // will be used as a source term in the
919 // adjoint eikonal equation
920 }
922 }
923
924 return tadjointEikonalSource;
925}
926
927
929{
930 const volVectorField& U = primalVars_.U();
931
933 const volTensorField& gradU = tgradU.cref();
934 tmp<volVectorField> tgradNuTilda = fvc::grad(nuTilda());
935 volVectorField& gradNuTilda = tgradNuTilda.ref();
936 tmp<volVectorField> tgradNuaTilda = fvc::grad(nuaTilda());
937 volVectorField::Boundary& gradNuaTildab =
938 tgradNuaTilda.ref().boundaryFieldRef();
939
940 // Explicitly correct the boundary gradient to get rid of
941 // the tangential component
942 forAll(mesh_.boundary(), patchI)
943 {
944 const fvPatch& patch = mesh_.boundary()[patchI];
945 if (isA<wallFvPatch>(patch))
946 {
947 tmp<vectorField> tnf = patch.nf();
948 const vectorField& nf = tnf();
949 // gradU:: can cause problems in zeroGradient patches for U
950 // and zero fixedValue for nuTilda.
951 // S becomes 0 and is used as a denominator in G
952 //gradU.boundaryField()[patchI] =
953 // nf * U_.boundaryField()[patchI].snGrad();
954 gradNuTilda.boundaryFieldRef()[patchI] =
955 nf*nuTilda().boundaryField()[patchI].snGrad();
956 gradNuaTildab[patchI] =
957 nf*nuaTilda().boundaryField()[patchI].snGrad();
958 }
959 }
960
961 // delta vorticity
962 volScalarField Omega(::sqrt(2.0)*mag(skew(gradU)));
963 volTensorField deltaOmega
964 (
965 (
966 (gradU & gradU)().T() //jk
967 - (gradU & gradU.T()) //symmetric
968 )
969 /(Omega + dimensionedScalar("SMALL", Omega.dimensions(), SMALL))
970 );
971 tgradU.clear();
972
973 volScalarField chi(this->chi());
974 volScalarField fv1(this->fv1(chi));
975 volScalarField fv2(this->fv2(chi, fv1));
976
981
982 // Assemply of the return field
983 auto tFISens
984 (
986 (
987 IOobject
988 (
989 type() + "flowTerm",
990 mesh_.time().timeName(),
991 mesh_,
994 ),
995 mesh_,
998 )
999 );
1000 volTensorField& FISens = tFISens.ref();
1001 FISens =
1002 // jk, cm formulation for the TM model convection
1003 - (nuaTilda()*(U*gradNuTilda))
1004 // jk, symmetric in theory
1005 + nuaTilda()*fvc::grad(DnuTildaEff()*gradNuTilda)().T()
1006 // jk
1007 - DnuTildaEff()*(tgradNuaTilda*gradNuTilda)
1008 // symmetric
1009 + 2.*nuaTilda()*Cb2_/sigmaNut_*(gradNuTilda*gradNuTilda)
1010 + (
1013 )
1014 *nuaTilda()*deltaOmega; // jk
1016
1017 return tFISens;
1018}
1019
1020
1022(
1023 const word& designVarsName
1024) const
1025{
1026 auto tres(tmp<scalarField>::New(nuTilda().primitiveField().size(), Zero));
1027 scalarField nuTildaSens
1028 (nuTilda().primitiveField()*nuaTilda().primitiveField());
1029
1032 (
1033 tres.ref(), nuTildaSens, fvOptions, nuTilda().name(), designVarsName
1034 );
1035
1036 return tres;
1037}
1038
1041{
1043}
1044
1045
1047{
1049 (
1051 "adjointSpalartAllmaras::correct"
1052 );
1053
1054 if (!adjointTurbulence_)
1055 {
1056 return;
1057 }
1058
1060
1062
1064 const volVectorField& Ua = adjointVars_.UaInst();
1065
1067 volScalarField gradUaR
1068 (
1070 );
1071
1072 dimensionedScalar oneOverSigmaNut = 1./sigmaNut_;
1073
1075
1077
1078 tmp<fvScalarMatrix> nuaTildaEqn
1079 (
1081 + fvm::div(-phi, nuaTilda())
1083 // Note: Susp
1085 + fvc::laplacian(2.0*Cb2_*oneOverSigmaNut*nuaTilda(), nuTilda())
1086 + gradNua*oneOverSigmaNut
1087 ==
1088 // Always a negative contribution to the lhs. No Sp used!
1090 // Always a positive contribution to the lhs. no need for SuSp
1092 - Cdnut_*gradUaR
1093 + fvOptions(nuaTilda())
1094 );
1095
1096 // Add sources from the objective functions
1097 objectiveManager_.addSource(nuaTildaEqn.ref());
1098
1099 nuaTildaEqn.ref().relax();
1100 fvOptions.constrain(nuaTildaEqn.ref());
1101 solve(nuaTildaEqn);
1103 nuaTilda().relax();
1104
1106 {
1107 dimensionedScalar maxDeltaNuaTilda =
1108 max(mag(nuaTilda() - nuaTilda().prevIter())());
1109 dimensionedScalar maxNuaTilda = max(mag(nuaTilda()));
1110 Info<< "Max mag (" << nuaTilda().name() << ") = "
1111 << maxNuaTilda.value() << endl;
1112 Info<< "Max mag (Delta " << nuaTilda().name() << ") = "
1113 << maxDeltaNuaTilda.value()
1114 << endl;
1115 }
1116}
1117
1118
1120{
1122 {
1125
1126 Cb1_.readIfPresent(this->coeffDict());
1127 Cb2_.readIfPresent(this->coeffDict());
1128 Cw1_ = Cb1_/sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_;
1129 Cw2_.readIfPresent(this->coeffDict());
1130 Cw3_.readIfPresent(this->coeffDict());
1131 Cv1_.readIfPresent(this->coeffDict());
1132 Cs_.readIfPresent(this->coeffDict());
1133
1134 return true;
1135 }
1136 else
1137 {
1138 return false;
1139 }
1140}
1141
1142
1143// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1144
1145} // End namespace adjointRASModels
1146} // End namespace incompressibleAdjoint
1147} // End namespace Foam
1148
1149// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const uniformDimensionedVectorField & g
fv::options & fvOptions
static tmp< volScalarField > createLimiter(const fvMesh &mesh, const dictionary &dict)
Return a limiter based on given cells. For use with classes other than ATC which employ the same smoo...
Definition ATCModel.C:216
label size() const noexcept
The number of elements in list.
Definition DLListBase.H:194
const dimensionSet & dimensions() const noexcept
Return dimensions.
void relax(const scalar alpha)
Relax field (for steady-state solution).
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())
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field).
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef().
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
void storePrevIter() const
Store the field as the previous iteration value.
void correctBoundaryConditions()
Correct boundary field.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
void setSize(label n)
Alias for resize().
Definition List.H:536
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
const word & name() const
Name function is needed to disambiguate those inherited from regIOobject and dictionary.
static autoPtr< boundaryAdjointContribution > New(const word &managerName, const word &adjointSolverName, const word &simulationType, const fvPatch &patch)
Return a reference to the selected turbulence model.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
Generic dimensioned Type class.
bool readIfPresent(const dictionary &dict)
Update the value of dimensioned<Type> if found in the dictionary, lookup in dictionary with the name(...
const Type & value() const noexcept
Return const reference to value.
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition fvMesh.H:395
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
Finite-volume options, which is an IOdictionary of values and a fv::optionList.
Definition fvOptions.H:69
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present otherwise lookup and return.
Definition fvOptions.C:116
Manages the adjoint mean flow fields and their mean values.
const volVectorField & Ua() const
Return const reference to velocity.
const volVectorField & UaInst() const
Return const reference to velocity.
const solverControl & getSolverControl() const
Return const reference to solverControl.
Abstract base class for incompressible turbulence models.
dictionary coeffDict_
Model coefficients dictionary.
bool includeDistance_
Does the turbulence model include distances and should the adjoint to the distance field be computed.
autoPtr< boundaryVectorField > wallShapeSensitivitiesPtr_
Wall sensitivity term for shape optimisation.
Switch adjointTurbulence_
Turbulence on/off flag.
bool changedPrimalSolution_
Has the primal solution changed?
autoPtr< boundaryVectorField > adjMomentumBCSourcePtr_
Source to the adjoint momentum BC emerging from differentiating the turbulence model.
const nearWallDist & y() const
Return the near wall distances.
autoPtr< volScalarField > adjointTMVariable1Ptr_
Adjoint turbulence model variable 1.
objectiveManager & objectiveManager_
Reference to the objectiveManager.
wordList adjointTMVariablesBaseNames_
Base names of the adjoint fields.
autoPtr< boundaryVectorField > wallFloCoSensitivitiesPtr_
Wall sensitivity term for flow control optimisation.
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
virtual bool read()
Read adjointRASProperties dictionary.
Continuous adjoint to the Spalart-Allmaras one-eqn mixing-length model for incompressible flows.
tmp< volScalarField > dfw_dNuTilda(const volScalarField &Stilda, const volScalarField &dfwdr, const volScalarField &dStildadNuTilda) const
virtual tmp< fvVectorMatrix > divDevReff(volVectorField &U) const
Return the diffusion term for the momentum equation.
virtual tmp< scalarField > diffusionCoeffVar1(label patchI) const
Diffusion coefficient of the first primal and adjoint turbulence model equation. Needed for some adjo...
const volScalarField & nuTilda() const
References to the primal turbulence model variables.
tmp< volScalarField > dStilda_dNuTilda(const volScalarField &Omega, const volScalarField &fv2, const volScalarField &dFv2dChi) const
virtual tmp< volTensorField > FISensitivityTerm()
Term contributing to the computation of FI-based sensitivities.
tmp< volScalarField > dfw_dOmega(const volScalarField &Stilda, const volScalarField &dfwdr, const volScalarField &dStildadOmega) const
tmp< volScalarField > dFv2_dChi(const volScalarField &chi, const volScalarField &fv1, const volScalarField &dFv1dChi) const
tmp< volScalarField > fw(const volScalarField &Stilda) const
tmp< volScalarField > dStilda_dDelta(const volScalarField &Omega, const volScalarField &fv2) const
virtual void correct()
Solve the adjoint turbulence equations.
tmp< volVectorField > conservativeMomentumSource()
Conservative source term for the adjoint momentum equations.
tmp< volScalarField > dStilda_dOmega(const volScalarField &Omega, const volScalarField &fv2) const
void updatePrimalRelatedFields()
Update the constant primal-related fields.
virtual const boundaryVectorField & adjointMomentumBCSource() const
Source for the outlet adjoint momentum BC coming from differentiating the turbulence model.
virtual tmp< volSymmTensorField > devReff() const
Return the effective stress tensor including the laminar stress.
tmp< volScalarField > dP_dNuTilda(const volScalarField &dStildadNuTilda) const
virtual void nullify()
Nullify all adjoint turbulence model fields and their old times.
tmp< volScalarField > dr_dNuTilda(const volScalarField &Stilda) const
tmp< volScalarField > fv1(const volScalarField &chi) const
virtual tmp< volScalarField > nutJacobianTMVar1() const
Jacobian of nut wrt the first turbulence model variable.
tmp< volScalarField > Stilda(const volScalarField &chi, const volScalarField &fv1) const
volScalarField & nuaTilda()
Access to the adjoint Spalart - Allmaras field.
tmp< volScalarField > dfw_dDelta(const volScalarField &Stilda, const volScalarField &dfwdr, const volScalarField &dStildadDelta) const
bool limitAdjointProduction_
Whether to limit the adjoint production term to enhance stability.
tmp< volScalarField > dfw_dr(const volScalarField &Stilda) const
virtual const boundaryVectorField & wallShapeSensitivities()
Sensitivity terms for shape optimisation, emerging from the turbulence model differentiation.
virtual tmp< volVectorField > adjointMeanFlowSource()
Source terms to the adjoint momentum equation due to the differentiation of the turbulence model.
virtual const boundaryVectorField & wallFloCoSensitivities()
Sensitivity terms for flow control, emerging from the turbulence model differentiation.
tmp< volScalarField > dD_dNuTilda(const volScalarField &fw, const volScalarField &dfwdNuTilda) const
virtual tmp< volScalarField > distanceSensitivities()
Sensitivity terms resulting from the differentiation of the distance field. Misses dxdb,...
tmp< volScalarField > dr_dStilda(const volScalarField &Stilda) const
tmp< volScalarField > dnut_dNuTilda(const volScalarField &fv1, const volScalarField &dFv1dChi) const
tmp< volScalarField > dr_dDelta(const volScalarField &Stilda) const
tmp< volScalarField > fv2(const volScalarField &chi, const volScalarField &fv1) const
volScalarField mask_
Field for masking (convergence enhancement).
virtual tmp< scalarField > topologySensitivities(const word &designVarsName) const
Term contributing to the computation of topology optimisation sensitivities.
tmp< volScalarField > dFv1_dChi(const volScalarField &chi) const
tmp< volScalarField > r(const volScalarField &Stilda) const
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
tmp< volScalarField > nu() const
Return the laminar viscosity.
virtual void correct()=0
Solve the adjoint turbulence equations.
Base class for solution control classes.
const volVectorField & U() const
Return const reference to velocity.
const surfaceScalarField & phi() const
Return const reference to volume flux.
This boundary condition provides a wall function for the turbulent viscosity (i.e....
Class for managing objective functions.
const word & adjointSolverName() const
Return name of the adjointSolver.
virtual void addSource(fvVectorMatrix &matrix)
Add contribution to adjoint momentum PDEs.
static void postProcessSens(scalarField &sens, scalarField &auxSens, fv::options &fvOptions, const word &fieldName, const word &designVariablesName)
Add part of the sensitivities coming from fvOptions.
bool printMaxMags() const
Print max mags of solver fields.
A class for managing temporary objects.
Definition tmp.H:75
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
Definition tmpI.H:221
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition tmpI.H:289
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
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 word & solverName() const
Return solver name.
bool useSolverNameForFields() const
Append solver name to fields?
static void nullifyField(GeometricField< Type, PatchField, GeoMesh > &fieldPtr)
Nullify field and old times, if present.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
const volScalarField & T
scalar uPlus
scalar yPlus
word timeName
Definition getTimeIndex.H:3
Namespace of functions to calculate explicit derivatives.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvmDiv.C:41
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvmDdt.C:41
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for incompressible adjoint turbulence models.
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedScalar pow6(const dimensionedScalar &ds)
dimensionedScalar pos(const dimensionedScalar &ds)
dimensionedScalar pow5(const dimensionedScalar &ds)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensionedScalar exp(const dimensionedScalar &ds)
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
GeometricField< tensor, fvPatchField, volMesh > volTensorField
dimensionedScalar sqrt(const dimensionedScalar &ds)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensionedScalar neg(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
volVectorField::Boundary boundaryVectorField
UList< label > labelUList
A UList of labels.
Definition UList.H:75
fvPatchField< vector > fvPatchVectorField
dimensionedTensor skew(const dimensionedTensor &dt)
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
fvPatchField< scalar > fvPatchScalarField
#define addProfiling(Name,...)
Define profiling trigger with specified name and description string. The description is generated by ...
volScalarField & nu
Calculation of adjoint based sensitivities for topology optimisation. This returns just the field par...
CEqn solve()
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299