adjointSpalartAllmaras.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) 2007-2019 PCOpt/NTUA
9 Copyright (C) 2013-2019 FOSS GP
10 Copyright (C) 2019 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
29Class
30 Foam::incompressibleAdjoint::adjointRASModels::adjointSpalartAllmaras
31
32Description
33 Continuous adjoint to the Spalart-Allmaras one-eqn mixing-length model for
34 incompressible flows.
35
36 Reference:
37 \verbatim
38 For the adjoint to the Spalart-Allmaras PDE
39
40 Zymaris, A., Papadimitriou, D., Giannakoglou, K., &
41 Othmer, C. (2009).
42 Continuous adjoint approach to the Spalart-Allmaras turbulence
43 model for incompressible flows.
44 Computers & Fluids, 38(8), 1528-538.
45 http://doi.org/10.1016/j.compfluid.2008.12.006
46
47 For the FI sensitivity terms
48
49 Papoutsis-Kiachagias, E. M., Asouti, V. G., Giannakoglou, K. C.,
50 Gkagkas, K., Shimokawa, S., & Itakura, E. (2018).
51 Multi-point aerodynamic shape optimization of cars based on
52 continuous adjoint.
53 Structural and Multidisciplinary Optimization, 59(2), 675-694.
54 http://doi.org/10.1007/s00158-018-2091-3
55
56 \endverbatim
57
58 Both of the above-mentioned papers develop the adjoint to the
59 Spalart-Allmaras PDE that includes the fv3 term. The current
60 implementation corresponds to the Spalart-Allmaras PDE as programmed within
61 OpenFOAM and is, thus, slightly different than the one developed in the
62 cited papers
63
64SourceFiles
65 adjointSpalartAllmaras.C
66
67\*---------------------------------------------------------------------------*/
68
69#ifndef adjointSpalartAllmaras_H
70#define adjointSpalartAllmaras_H
71
72#include "adjointRASModel.H"
73#include "wallDist.H"
74
75// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
76
77namespace Foam
78{
79namespace incompressibleAdjoint
80{
81namespace adjointRASModels
82{
83
84/*---------------------------------------------------------------------------*\
85 Class adjointSpalartAllmaras Declaration
86\*---------------------------------------------------------------------------*/
89:
90 public adjointRASModel
91{
92 // Private Member Functions
93
94 //- No copy construct
96
97 //- No copy assignment
98 void operator=(const adjointSpalartAllmaras&) = delete;
99
100
101protected:
102
103 // Protected data
104
105 // Model coefficients
117
118 //- Whether to limit the adjoint production term to enhance
119 //- stability
121
122
123 // Fields
124
125 //- Wall distance
126 // Note: reference to the distance known by the primal model
127 const volScalarField& y_;
128
129 //- Field for masking (convergence enhancement)
131
132 // Fields that depend only on primal fields and are very expensive
133 // to compute in each iteration.
134 // Store and update them when the primal solution has been updated
144
145 // Useful quantities for bounding
147
148
149 // Protected Member Functions
150
151 // Primal Spalart - Allmaras
152
153 tmp<volScalarField> chi() const;
154
156
158 (
159 const volScalarField& chi,
160 const volScalarField& fv1
161 ) const;
162
164 (
165 const volScalarField& chi,
166 const volScalarField& fv1
167 ) const;
168
170
172
174
175 //- References to the primal turbulence model variables
176 const volScalarField& nuTilda() const;
177
178 const volScalarField& nut() const;
179
180
181 // Adjoint Spalart - Allmaras
182
183 // Differentiation of the primal Spalart - Allmaras terms
184
186 (
187 const volScalarField& chi
188 ) const;
189
191 (
192 const volScalarField& chi,
193 const volScalarField& fv1,
194 const volScalarField& dFv1dChi
195 ) const;
196
198 (
199 const volScalarField& Omega,
200 const volScalarField& fv2
201 ) const;
202
204 (
205 const volScalarField& Omega,
206 const volScalarField& fv2,
207 const volScalarField& dFv2dChi
208 ) const;
209
211 (
212 const volScalarField& Omega,
213 const volScalarField& fv2
214 ) const;
215
217 (
219 ) const;
220
222 (
224 ) const;
225
227 (
229 ) const;
230
232 (
234 ) const;
235
237 (
238 const volScalarField& Stilda,
239 const volScalarField& dfwdr,
240 const volScalarField& dStildadNuTilda
241 ) const;
242
244 (
245 const volScalarField& Stilda,
246 const volScalarField& dfwdr,
247 const volScalarField& dStildadOmega
248 ) const;
249
251 (
252 const volScalarField& Stilda,
253 const volScalarField& dfwdr,
254 const volScalarField& dStildadDelta
255 ) const;
256
258 (
259 const volScalarField& fw,
260 const volScalarField& dfwdNuTilda
261 ) const;
262
264 (
265 const volScalarField& dStildadNuTilda
266 ) const;
267
269 (
270 const volScalarField& fv1,
271 const volScalarField& dFv1dChi
272 ) const;
273
274
275 //- Conservative source term for the adjoint momentum equations
276 // Sets also the adjointMomentumBCSource
278
279 //- Access to the adjoint Spalart - Allmaras field
280 inline volScalarField& nuaTilda()
281 {
282 return adjointTMVariable1Ptr_();
283 };
284
285 //- Update the constant primal-related fields
287
288 //- Allocate the mask field
290
291
292public:
293
294 //- Runtime type information
295 TypeName("adjointSpalartAllmaras");
296
297
298 // Constructors
299
300 //- Construct from components
302 (
303 incompressibleVars& primalVars,
305 objectiveManager& objManager,
306 const word& adjointTurbulenceModelName
308 const word& modelName = typeName
309 );
310
311
312 //- Destructor
313 virtual ~adjointSpalartAllmaras() = default;
314
315
316 // Member Functions
317
318 virtual tmp<volSymmTensorField> devReff() const;
319
320 virtual tmp<volSymmTensorField> devReff(const volVectorField& U) const;
321
323
325
327
328 virtual tmp<scalarField> diffusionCoeffVar1(label patchI) const;
329
330 virtual const boundaryVectorField& adjointMomentumBCSource() const;
331
333
335
337
339
340 //- Nullify all adjoint turbulence model fields and their old times
341 virtual void nullify();
342
343 //- Solve the adjoint turbulence equations
344 virtual void correct();
345
346 //- Read adjointRASProperties dictionary
347 virtual bool read();
348
349};
350
351
352// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
353
354} // End namespace adjointRASModels
355} // End namespace incompressibleAdjoint
356} // End namespace Foam
357
358// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
359
360#endif
361
362// ************************************************************************* //
Generic GeometricBoundaryField class.
Manages the adjoint mean flow fields and their mean values.
Abstract base class for incompressible turbulence models.
autoPtr< volScalarField > adjointTMVariable1Ptr_
Adjoint turbulence model variable 1.
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.
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
TypeName("adjointSpalartAllmaras")
Runtime type information.
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 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
tmp< volScalarField > dfw_dr(const volScalarField &Stilda) const
virtual const boundaryVectorField & wallShapeSensitivities()
Sensitivity terms for shape optimisation, emerging from.
virtual const boundaryVectorField & wallFloCoSensitivities()
Sensitivity terms for flow control, emerging from the.
tmp< volScalarField > dD_dNuTilda(const volScalarField &fw, const volScalarField &dfwdNuTilda) const
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
const volScalarField & nut() const
Return the turbulence viscosity.
tmp< volScalarField > fv2(const volScalarField &chi, const volScalarField &fv1) const
volScalarField mask_
Field for masking (convergence enhancement)
tmp< volScalarField > dFv1_dChi(const volScalarField &chi) const
tmp< volScalarField > r(const volScalarField &Stilda) const
Base class for solution control classes.
class for managing incompressible objective functions.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
Namespace for OpenFOAM.
static const char *const typeName
The type name used in ensight case files.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73