Loading...
Searching...
No Matches
objectiveNutSqr.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-2023 PCOpt/NTUA
9 Copyright (C) 2013-2023 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::objectives::objectiveNutSqr
30
31Description
32 Objective qualitatively quantifying noise through the integral of the
33 squared turbulent viscosity over specified cellZones. Requires the adjoint
34 to the turbulence model to be incorporated into the optimisation loop.
35
36 Reference:
37 \verbatim
38 Objective function initially presented in
39
40 Papoutsis-Kiachagias, E. M., Magoulas, N., Mueller, J., Othmer, C.,
41 & Giannakoglou, K. C. (2015).
42 Noise reduction in car aerodynamics using a surrogate objective
43 function and the continuous adjoint method with wall functions.
44 Computers & Fluids, 122(20), 223-232.
45 https://doi.org/10.1016/j.compfluid.2015.09.002
46 \endverbatim
47
48SourceFiles
49 objectiveNutSqrNoise.C
50
51\*---------------------------------------------------------------------------*/
52
53#ifndef objectiveNutSqr_H
54#define objectiveNutSqr_H
55
57
58// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59
60namespace Foam
61{
62
63namespace objectives
64{
66/*---------------------------------------------------------------------------*\
67 Class objectiveNutSqr Declaration
68\*---------------------------------------------------------------------------*/
69
71:
72 public objectiveIncompressible
73{
74 // Private data
75
76 //- Where to define the objective
77 labelList zones_;
78
79 //- List with the names of the adjoint turbulence model fields
80 // This is kept separately from fieldNames since the latter may
81 // or may not include Ua, depending on whether nut is a function
82 // of U. This makes deciding on whether to add or not sources
83 // to a given fvScalarMatrix tricky, hence the utilisation of
84 // adjointTurbulenceNames_
85 wordList adjointTurbulenceNames_;
86
87
88 // Private Member Functions
89
90 //- Populate fieldNames
91 void populateFieldNames();
92
93
94public:
95
96 //- Runtime type information
97 TypeName("nutSqr");
98
99
100 // Constructors
102 //- Construct from components
104 (
105 const fvMesh& mesh,
106 const dictionary& dict,
107 const word& adjointSolverName,
108 const word& primalSolverName
109 );
110
111
112 //- Destructor
113 virtual ~objectiveNutSqr() = default;
114
115
116 // Member Functions
117
118 //- Return the objective function value
119 virtual scalar J();
120
121 //- Update values to be added to the adjoint outlet velocity
122 virtual void update_dJdv();
123
124 //- Update field to be added to the first adjoint turbulence model PDE
125 virtual void update_dJdTMvar1();
126
127 //- Update field to be added to the second adjoint turbulence model PDE
128 virtual void update_dJdTMvar2();
129
130 //- Update field to be added to be added to volume-based
131 //- sensitivity derivatives, emerging from delta ( dV ) / delta b
132 virtual void update_divDxDbMultiplier();
133
134 //- Add source terms to the adjoint turbulence model equations
135 virtual void addSource(fvScalarMatrix& matrix);
136};
137
138
139// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
140
141} // End namespace objectives
142} // End namespace Foam
143
144// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
145
146#endif
147
148// ************************************************************************* //
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const dictionary & dict() const
Return objective dictionary.
Definition objective.C:91
virtual void update_dJdTMvar2()
Update field to be added to the second adjoint turbulence model PDE.
objectiveNutSqr(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
Construct from components.
virtual void addSource(fvScalarMatrix &matrix)
Add source terms to the adjoint turbulence model equations.
virtual scalar J()
Return the objective function value.
virtual void update_dJdTMvar1()
Update field to be added to the first adjoint turbulence model PDE.
TypeName("nutSqr")
Runtime type information.
virtual void update_divDxDbMultiplier()
Update field to be added to be added to volume-based sensitivity derivatives, emerging from delta ( d...
virtual ~objectiveNutSqr()=default
Destructor.
virtual void update_dJdv()
Update values to be added to the adjoint outlet velocity.
A class for handling words, derived from Foam::string.
Definition word.H:66
dynamicFvMesh & mesh
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
List< label > labelList
A List of labels.
Definition List.H:62
fvMatrix< scalar > fvScalarMatrix
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68