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-2022 PCOpt/NTUA
9  Copyright (C) 2013-2022 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 Class
29  Foam::objectives::objectiveNutSqr
30 
31 Description
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 
48 SourceFiles
49  objectiveNutSqrNoise.C
50 
51 \*---------------------------------------------------------------------------*/
52 
53 #ifndef objectiveNutSqr_H
54 #define objectiveNutSqr_H
55 
57 
58 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59 
60 namespace Foam
61 {
62 
63 namespace objectives
64 {
65 
66 /*---------------------------------------------------------------------------*\
67  Class objectiveNutSqr Declaration
68 \*---------------------------------------------------------------------------*/
69 
70 class objectiveNutSqr
71 :
73 {
74  // Private data
75 
76  //- Where to define the objective
77  labelList zones_;
78 
79 
80 public:
81 
82  //- Runtime type information
83  TypeName("nutSqr");
84 
85 
86  // Constructors
87 
88  //- Construct from components
90  (
91  const fvMesh& mesh,
92  const dictionary& dict,
93  const word& adjointSolverName,
94  const word& primalSolverName
95  );
96 
97 
98  //- Destructor
99  virtual ~objectiveNutSqr() = default;
100 
101 
102  // Member Functions
103 
104  //- Return the objective function value
105  scalar J();
106 
107  //- Update values to be added to the adjoint outlet velocity
108  void update_dJdv();
109 
110  //- Update field to be added to the first adjoint turbulence model PDE
111  void update_dJdTMvar1();
112 
113  //- Update field to be added to the second adjoint turbulence model PDE
114  void update_dJdTMvar2();
115 
116  //- Update field to be added to be added to volume-based
117  //- sensitivity derivatives, emerging from delta ( dV ) / delta b
119 };
120 
121 
122 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
123 
124 } // End namespace objectives
125 } // End namespace Foam
126 
127 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
128 
129 #endif
130 
131 // ************************************************************************* //
objectiveNutSqr(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
Construct from components.
virtual ~objectiveNutSqr()=default
Destructor.
const dictionary & dict() const
Return objective dictionary.
Definition: objective.C:87
void update_dJdTMvar2()
Update field to be added to the second adjoint turbulence model PDE.
Objective qualitatively quantifying noise through the integral of the squared turbulent viscosity ove...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
void update_dJdv()
Update values to be added to the adjoint outlet velocity.
TypeName("nutSqr")
Runtime type information.
scalar J()
Return the objective function value.
void update_divDxDbMultiplier()
Update field to be added to be added to volume-based sensitivity derivatives, emerging from delta ( d...
dynamicFvMesh & mesh
A class for handling words, derived from Foam::string.
Definition: word.H:63
void update_dJdTMvar1()
Update field to be added to the first adjoint turbulence model PDE.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
Abstract base class for objective functions in incompressible flows.
Namespace for OpenFOAM.