Loading...
Searching...
No Matches
adjointkOmegaSST.H
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-2022 PCOpt/NTUA
9 Copyright (C) 2014-2022 FOSS GP
10-------------------------------------------------------------------------------
11License
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
28Class
29 Foam::incompressibleAdjoint::adjointRASModels::adjointkOmegaSST
30
31Description
32 Continuous adjoint to the kOmegaSST turbulence model for incompressible
33 flows.
34
35 Reference:
36 \verbatim
37 The code is based on the following reference, with a number of
38 changes in the numerical implementation
39
40 Kavvadias, I., Papoutsis-Kiachagias, E., Dimitrakopoulos, G., &
41 Giannakoglou, K. (2014).
42 The continuous adjoint approach to the k–ω SST turbulence model
43 with applications in shape optimization
44 Engineering Optimization, 47(11), 1523-1542.
45 https://doi.org/10.1080/0305215X.2014.979816
46 \endverbatim
47
48SourceFiles
49 adjointkOmegaSST.C
50
51\*---------------------------------------------------------------------------*/
52
53#ifndef Foam_adjointkOmegaSST_H
54#define Foam_adjointkOmegaSST_H
55
56#include "adjointRASModel.H"
57#include "wallDist.H"
58
59// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60
61namespace Foam
62{
64{
65namespace adjointRASModels
66{
68/*---------------------------------------------------------------------------*\
69 Class adjointkOmegaSST Declaration
70\*---------------------------------------------------------------------------*/
71
72class adjointkOmegaSST
73:
74 public adjointRASModel
75{
76 // Private Member Functions
77
78 //- No copy construct
79 adjointkOmegaSST(const adjointkOmegaSST&) = delete;
80
81 //- No copy assignment
82 void operator=(const adjointkOmegaSST&) = delete;
83
84
85protected:
86
87 // Protected data
88
89 // Primal Model coefficients
109
110 //- Flag to include the F3 term
111 Switch F3_;
113
114 // Fields
115
116 // Primal Fields
117
118 //- Wall distance
119 // Note: reference to the distance known by the primal model
120 const volScalarField& y_;
121
122 //- Cached primal gradient fields
126
127 //- Primal cached fields involved in the solution of the
128 // adjoint equations
129 // Cached to reduce the computational cost
137
138 // Model Field coefficients
145 // Switch fields
146
147 // Switch fields for the differentiation of F1
152
153 //- Switch fields for the production in the k Eqn
158 // Switch fields for the differentiation of nut
159 // Holds also for the differentiation of the second branch of
160 // GbyNu
164
165 // Switch fields for GPrime
168
169 // Zero first cell field and IDs
170 // Since the omega equation is a two-zonal one, some
171 // of the terms in the adjoint equations need to ba canceled
172 // at the cells next to omegaWallFunction BCs
175
176
177 // Turbulence model multipliers
179 //- Nut Jacobian w.r.t. omega
181
182 //- Nut Jacobian w.r.t. k
185 //- Diffusivity of the omega equation
187
188 //- Diffusivity of the k equation
190
191
192 // Protected Member Functions
194 // Primal functions
195
196 virtual tmp<volScalarField> F1() const;
197 virtual tmp<volScalarField> F2() const;
199 (
200 const volScalarField& GbyNu0,
201 const volScalarField& F2,
202 const volScalarField& S2
203 ) const;
204
205 //- Return G/nu
207 (
211 ) const;
212
214 (
215 const volScalarField& F1,
216 const dimensionedScalar& psi1,
218 ) const
219 {
220 return F1*(psi1 - psi2) + psi2;
221 }
222
224 (
226 const dimensionedScalar& psi1,
228 ) const
229 {
230 return F1*(psi1 - psi2) + psi2;
231 }
232
235 return blend(F1, alphaK1_, alphaK2_);
236 }
237
239 {
241 }
242
246 ) const
247 {
249 (
250 IOobject::scopedName(this->type(), "beta"),
251 blend(F1, beta1_, beta2_)
252 );
253 }
258 ) const
261 (
262 IOobject::scopedName(this->type(), "beta"),
263 blend(F1, beta1_, beta2_)
264 );
265 }
266
268 (
270 ) const
271 {
273 (
274 IOobject::scopedName(this->type(), "gamma"),
275 blend(F1, gamma1_, gamma2_)
276 );
277 }
278
280 (
281 const volScalarField& F1
282 ) const
283 {
285 (
286 IOobject::scopedName(this->type(), "gamma"),
287 blend(F1, gamma1_, gamma2_)
288 );
289 }
290
292
293
294 // References to the primal turbulence model variables
295
296 inline const volScalarField& k() const
297 {
298 return primalVars_.RASModelVariables()().TMVar1();
299 }
301 inline volScalarField& k()
302 {
303 return primalVars_.RASModelVariables()().TMVar1();
304 }
305
306 inline const volScalarField& omega() const
307 {
308 return primalVars_.RASModelVariables()().TMVar2();
309 }
310
311 inline volScalarField& omega()
312 {
313 return primalVars_.RASModelVariables()().TMVar2();
314 }
315
316 inline const volScalarField& nutRef() const
318 return primalVars_.RASModelVariables()().nutRef();
319 }
320
321 inline volScalarField& nutRef()
323 return primalVars_.RASModelVariables()().nutRef();
324 }
325
326
327 // Adjoint related protected member functions
328
329 //- Derivative of the primal equations wrt nut
331
332 //- Nut Jacobian wrt omega
334 (
335 const volScalarField& F2,
336 const volScalarField& S,
337 const volScalarField& case_1_nut,
338 const volScalarField& case_2_nut,
339 const volScalarField& case_3_nut
340 ) const;
341
342 //- Nut Jacobian wrt k
344 (
345 const volScalarField& F2,
346 const volScalarField& S,
347 const volScalarField& case_2_nut
348 ) const;
349
350 //- F2 Jacobian wrt omega
352 (
353 const volScalarField& F2,
354 const volScalarField& case_2_nut,
355 const volScalarField& case_3_nut
356 ) const;
357
358 //- F2 Jacobian wrt k
360 (
361 const volScalarField& F2,
362 const volScalarField& case_2_nut
363 ) const;
364
365 //- GbyNu Jacobian wrt omega
367
368 //- GbyNu Jacobian wrt k
370
371 //- Derivative of the primal equations wrt F1
373
374 //- F1 Jacobian wrt omega (no contributions from grad(omega))
376
377 //- F1 Jacobian wrt grad(omega)
379 (
380 const volScalarField& arg1
381 ) const;
382
383 //- Source to waEqn from the differentiation of F1
385
386 //- Source to waEqn from the differentiation of CDkOmega
388
389 //- F1 Jacobian wrt k (no contributions from grad(k))
390 tmp<volScalarField> dF1_dk(const volScalarField& arg1) const;
391
392 //- F1 Jacobian wrt grad(k)
394
395 //- Source to kaEqn from the differentiation of F1
397
398 //- Source to kaEqn from the differentiation of CDkOmega
400
401 //- Differentiation of the turbulence model diffusion coefficients
403 (
404 const volScalarField& primalField,
405 const volScalarField& adjointField,
406 const word& schemeName
407 ) const;
408
409 //- Term multiplying dnut/db, coming from the turbulence model
411 (
412 const volScalarField& primalField,
413 const volScalarField& adjointField,
414 const volScalarField& coeffField,
415 const volScalarField& bcField,
416 const word& schemeName
417 ) const;
418
419 //- Term multiplying dnut/db, coming from the momentum equations
421 (
422 const volVectorField& primalField,
423 const volVectorField& adjointField,
424 const volScalarField& bcField,
425 const word& schemeName
426 ) const;
427
428 // Functions computing the adjoint mean flow source
429
430 //- Contributions from the turbulence model convection terms
432 (
433 const volScalarField& primalField,
434 const volScalarField& adjointField
435 ) const;
436
437 //- Contributions from the G
439 (
440 tmp<volSymmTensorField>& GbyNuMult
441 ) const;
442
443 //- Contributions from the divU
445 (
446 tmp<volScalarField>& divUMult
447 ) const;
448
449 //- Contributions from nut(U), in the diffusion coefficients
450 //- of the turbulence model
452 (
453 const volScalarField& primalField,
454 const volScalarField& adjointField,
455 const volScalarField& coeffField
456 ) const;
457
458 //- Contributions from nut(U)
460 (
462 const volScalarField& F2,
463 const volScalarField& S,
464 const volScalarField& case_1_nut,
465 const volTensorField& gradU
466 ) const;
467
468
469 //- Contributions from the differentiation of k existing in
470 //- nutkWallFunction.
471 // This could also be implemented in kaqRWallFunction but all
472 // the fields required for the computation already exist here,
473 // hence the code complexity is reduced
475 (
476 fvScalarMatrix& kaEqn,
478 );
479
480 // References to the adjoint turbulence model fields
481
482 inline volScalarField& ka()
483 {
484 return adjointTMVariable1Ptr_();
485 };
486
487 inline const volScalarField& ka() const
488 {
489 return adjointTMVariable1Ptr_();
490 };
491
492 inline volScalarField& wa()
493 {
494 return adjointTMVariable2Ptr_();
495 };
496
497 inline const volScalarField& wa() const
498 {
499 return adjointTMVariable2Ptr_();
500 };
501
502
503 //- Update of the primal cached fields
505
506 //- Return the requested interpolation scheme if it exists,
507 //- otherwise return a reverseLinear scheme
508 template<class Type>
509 tmp<surfaceInterpolationScheme<Type>> interpolationScheme
510 (
511 const word& schemeName
512 ) const;
513
514 //- Return the interpolation scheme used by the primal convection
515 //- term of the equation corresponding to the argument
516 tmp<surfaceInterpolationScheme<scalar>> convectionScheme
517 (
518 const word& varName
519 ) const;
520
521
522public:
523
524 //- Runtime type information
525 TypeName("adjointkOmegaSST");
526
527
528 // Constructors
529
530 //- Construct from components
531 adjointkOmegaSST
532 (
533 incompressibleVars& primalVars,
534 incompressibleAdjointMeanFlowVars& adjointVars,
535 objectiveManager& objManager,
536 const word& adjointTurbulenceModelName
537 = adjointTurbulenceModel::typeName,
538 const word& modelName = typeName
539 );
540
541
542 //- Destructor
543 virtual ~adjointkOmegaSST() = default;
544
545
546 // Member Functions
547
548 //- Return the effective diffusivity for k
549 tmp<volScalarField> DkEff(const volScalarField& F1) const
550 {
551 return tmp<volScalarField>
552 (
553 new volScalarField("DkEff", alphaK(F1)*this->nut() + this->nu())
554 );
555 }
556
557 //- Return the effective diffusivity for omega
561 (
563 (
564 "DomegaEff",
565 alphaOmega(F1)*this->nut() + this->nu()
566 )
567 );
568 }
570 virtual tmp<volSymmTensorField> devReff() const;
571
572 virtual tmp<volSymmTensorField> devReff(const volVectorField& U) const;
573
574 //- Return the transpose part of the adjoint momentum stresses
576
577 //- Non-conservative part of the terms added to the mean flow equations
579
580 //- Source term added to the adjoint mean flow due to the
581 //- differentiation of the turbulence model
583
584 //- Jacobian of nut wrt to k
585 // Needs to be implemented for objectives related to nut, defined in
586 // the internal field
588
589 //- Jacobian of nut wrt to omega
590 // Needs to be implemented for objectives related to nut, defined in
591 // the internal field
593
594 //- Jacobian of nut wrt the flow velocity
595 // Assumes we want to get contributions of mult*d(nut)/dU
596 // Since the dependency of nut to U is usually through a differential
597 // operator, the term multiplying d(nut)/dU is passed as an argument
598 // to this function; the latter should then compute the
599 // contribution of the afforementioned term to the adjoint mean flow
600 // equations
602 (
603 tmp<volScalarField>& dNutdUMult
604 ) const;
606 //- Diffusion coeff at the boundary for k
607 virtual tmp<scalarField> diffusionCoeffVar1(label patchI) const;
608
609 //- Diffusion coeff at the boundary for omega
610 virtual tmp<scalarField> diffusionCoeffVar2(label patchI) const;
611
612 virtual const boundaryVectorField& adjointMomentumBCSource() const;
613
614 //- Sensitivity derivative contributions when using the (E)SI approach
616
618
619 //- Contributions to the adjoint eikonal equation (zero for now)
621
622 //- Sensitivity derivative contributions when using the FI approach
624
626 (
627 const word& designVarsName
628 ) const;
629
630 //- Nullify all adjoint turbulence model fields and their old times
631 virtual void nullify();
632
633 //- Solve the adjoint turbulence equations
634 virtual void correct();
636 //- Read adjointRASProperties dictionary
637 virtual bool read();
638};
639
640
641// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
642
643} // End namespace adjointRASModels
644} // End namespace incompressibleAdjoint
645} // End namespace Foam
647// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
648
649#ifdef NoRepository
651#endif
652
653// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
654
655#endif
656
657// ************************************************************************* //
#define F2(B, C, D)
Definition SHA1.C:150
#define F1(B, C, D)
Definition SHA1.C:149
const volScalarField & psi2
const volScalarField & psi1
DimensionedField< scalar, volMesh > Internal
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition Switch.H:81
Manages the adjoint mean flow fields and their mean values.
autoPtr< volScalarField > adjointTMVariable1Ptr_
Adjoint turbulence model variable 1.
autoPtr< volScalarField > adjointTMVariable2Ptr_
Adjoint turbulence model variable 2.
virtual tmp< fvVectorMatrix > divDevReff(volVectorField &U) const
Return the transpose part of the adjoint momentum stresses.
virtual tmp< scalarField > diffusionCoeffVar1(label patchI) const
Diffusion coeff at the boundary for k.
virtual tmp< volScalarField > nutJacobianTMVar2() const
Jacobian of nut wrt to omega.
tmp< volScalarField > dGPrime_domega() const
GbyNu Jacobian wrt omega.
tmp< volVectorField > GMeanFlowSource(tmp< volSymmTensorField > &GbyNuMult) const
Contributions from the G.
virtual tmp< volTensorField > FISensitivityTerm()
Sensitivity derivative contributions when using the FI approach.
tmp< volVectorField > dF1_dGradK(const volScalarField &arg1) const
F1 Jacobian wrt grad(k).
tmp< volScalarField::Internal > blend(const volScalarField::Internal &F1, const dimensionedScalar &psi1, const dimensionedScalar &psi2) const
tmp< volScalarField > dR_dF1() const
Derivative of the primal equations wrt F1.
tmp< volScalarField > dF1_dk(const volScalarField &arg1) const
F1 Jacobian wrt k (no contributions from grad(k)).
tmp< volScalarField > dGPrime_dk() const
GbyNu Jacobian wrt k.
tmp< volVectorField > convectionMeanFlowSource(const volScalarField &primalField, const volScalarField &adjointField) const
Contributions from the turbulence model convection terms.
tmp< volScalarField::Internal > gamma(const volScalarField::Internal &F1) const
tmp< surfaceInterpolationScheme< Type > > interpolationScheme(const word &schemeName) const
Return the requested interpolation scheme if it exists, otherwise return a reverseLinear scheme.
tmp< volScalarField > dF2_dk(const volScalarField &F2, const volScalarField &case_2_nut) const
F2 Jacobian wrt k.
virtual void correct()
Solve the adjoint turbulence equations.
volScalarField DkEff_
Diffusivity of the k equation.
volScalarField S2_
Primal cached fields involved in the solution of the.
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.
tmp< fvScalarMatrix > waEqnSourceFromCDkOmega() const
Source to waEqn from the differentiation of CDkOmega.
tmp< volScalarField > kaEqnSourceFromF1() const
Source to kaEqn from the differentiation of F1.
tmp< volScalarField > dR_dnut()
Derivative of the primal equations wrt nut.
void updatePrimalRelatedFields()
Update of the primal cached fields.
virtual const boundaryVectorField & adjointMomentumBCSource() const
Source for the outlet adjoint momentum BC coming from differentiating the turbulence model.
tmp< surfaceInterpolationScheme< scalar > > convectionScheme(const word &varName) const
Return the interpolation scheme used by the primal convection term of the equation corresponding to t...
virtual tmp< volSymmTensorField > devReff() const
Return the effective stress tensor including the laminar stress.
tmp< volScalarField > beta(const volScalarField &F1) const
virtual void nullify()
Nullify all adjoint turbulence model fields and their old times.
tmp< volScalarField > blend(const volScalarField &F1, const dimensionedScalar &psi1, const dimensionedScalar &psi2) const
virtual tmp< volScalarField > nutJacobianTMVar1() const
Jacobian of nut wrt to k.
tmp< volScalarField > dnut_dk(const volScalarField &F2, const volScalarField &S, const volScalarField &case_2_nut) const
Nut Jacobian wrt k.
volTensorField gradU_
Cached primal gradient fields.
volScalarField DOmegaEff_
Diffusivity of the omega equation.
virtual tmp< volScalarField > GbyNu(const volScalarField &GbyNu0, const volScalarField &F2, const volScalarField &S2) const
tmp< volScalarField > alphaK(const volScalarField &F1) const
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.
virtual const boundaryVectorField & wallShapeSensitivities()
Sensitivity derivative contributions when using the (E)SI approach.
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).
tmp< volScalarField > DomegaEff(const volScalarField &F1) const
Return the effective diffusivity for omega.
tmp< volScalarField > alphaOmega(const volScalarField &F1) const
void addWallFunctionTerms(fvScalarMatrix &kaEqn, const volScalarField &dR_dnut)
Contributions from the differentiation of k existing in nutkWallFunction.
tmp< volScalarField > waEqnSourceFromF1() const
Source to waEqn from the differentiation of F1.
virtual tmp< volVectorField > adjointMeanFlowSource()
Source term added to the adjoint mean flow due to the differentiation of the turbulence model.
volScalarField case_1_Pk_
Switch fields for the production in the k Eqn.
tmp< volScalarField::Internal > beta(const volScalarField::Internal &F1) const
virtual const boundaryVectorField & wallFloCoSensitivities()
Sensitivity terms for flow control, emerging from the turbulence model differentiation.
virtual tmp< volVectorField > nonConservativeMomentumSource() const
Non-conservative part of the terms added to the mean flow equations.
TypeName("adjointkOmegaSST")
Runtime type information.
tmp< volScalarField > kaEqnSourceFromCDkOmega() const
Source to kaEqn from the differentiation of CDkOmega.
tmp< volScalarField > coeffsDifferentiation(const volScalarField &primalField, const volScalarField &adjointField, const word &schemeName) const
Differentiation of the turbulence model diffusion coefficients.
virtual tmp< volScalarField > distanceSensitivities()
Contributions to the adjoint eikonal equation (zero for now).
virtual tmp< volVectorField > nutJacobianU(tmp< volScalarField > &dNutdUMult) const
Jacobian of nut wrt the flow velocity.
virtual tmp< scalarField > diffusionCoeffVar2(label patchI) const
Diffusion coeff at the boundary for omega.
tmp< volVectorField > divUMeanFlowSource(tmp< volScalarField > &divUMult) const
Contributions from the divU.
tmp< volVectorField > dF1_dGradOmega(const volScalarField &arg1) const
F1 Jacobian wrt grad(omega).
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.
virtual tmp< scalarField > topologySensitivities(const word &designVarsName) const
Term contributing to the computation of topology optimisation sensitivities.
tmp< volScalarField > DkEff(const volScalarField &F1) const
Return the effective diffusivity for k.
tmp< volScalarField > dF1_domega(const volScalarField &arg1) const
F1 Jacobian wrt omega (no contributions from grad(omega)).
virtual bool read()
Read adjointRASProperties dictionary.
tmp< volScalarField > gamma(const volScalarField &F1) const
tmp< volScalarField > dF2_domega(const volScalarField &F2, const volScalarField &case_2_nut, const volScalarField &case_3_nut) const
F2 Jacobian wrt omega.
tmp< volScalarField > nu() const
Return the laminar viscosity.
virtual const tmp< volScalarField > nut() const
Return the turbulence viscosity.
Base class for solution control classes.
const autoPtr< incompressible::RASModelVariables > & RASModelVariables() const
Return const reference to the turbulence model variables.
Class for managing objective functions.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
const scalar gamma
Definition EEqn.H:9
Namespace for incompressible adjoint turbulence models.
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
List< label > labelList
A List of labels.
Definition List.H:62
fvMatrix< scalar > fvScalarMatrix
GeometricField< scalar, fvPatchField, volMesh > volScalarField
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
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
GeometricField< tensor, fvPatchField, volMesh > volTensorField
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
volVectorField::Boundary boundaryVectorField
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68