objectivePowerDissipation.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::objectivePowerDissipation
30 
31 Description
32  Computes and minimizes the power dissipation within given cellZones.
33  In the absence of significant viscous stresses on the "inlet" and "outlet"
34  of the cellZones, this value is equal to the volume flow rate-weigthed
35  total pressure losses (see also objectivePtLosses) within th cellZones
36 
37 SourceFiles
38  objectivePowerDissipation.C
39 
40 \*---------------------------------------------------------------------------*/
41 
42 #ifndef objectivePowerDissipation_H
43 #define objectivePowerDissipation_H
44 
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 namespace objectives
53 {
54 
55 /*---------------------------------------------------------------------------*\
56  Class objectivePowerDissipation Declaration
57 \*---------------------------------------------------------------------------*/
58 
60 :
62 {
63  // Private data
64 
65  //- Where to define the objective
66  labelList zones_;
67 
68 
69 public:
70 
71  //- Runtime type information
72  TypeName("powerDissipation");
73 
74 
75  // Constructors
76 
77  //- Construct from components
79  (
80  const fvMesh& mesh,
81  const dictionary& dict,
82  const word& adjointSolverName,
83  const word& primalSolverName
84  );
85 
86 
87  //- Destructor
88  virtual ~objectivePowerDissipation() = default;
89 
90 
91  // Member Functions
92 
93  //- Return the objective function value
94  scalar J();
95 
96  //- Update values to be added to the adjoint outlet velocity
97  void update_dJdv();
98 
99  //- Update field to be added to the first adjoint turbulence model PDE
100  void update_dJdTMvar1();
101 
102  //- Update field to be added to the second adjoint turbulence model PDE
103  void update_dJdTMvar2();
104 
105  //- Update div(dx/db multiplier). Volume-based sensitivity term
106  virtual void update_divDxDbMultiplier();
107 
108  //- Update grad(dx/db multiplier). Volume-based sensitivity term
109  virtual void update_gradDxDbMultiplier();
110 };
111 
112 
113 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
114 
115 } // End namespace objectives
116 } // End namespace Foam
117 
118 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
119 
120 #endif
121 
122 // ************************************************************************* //
objectivePowerDissipation(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
Construct from components.
const dictionary & dict() const
Return objective dictionary.
Definition: objective.C:87
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.
scalar J()
Return the objective function value.
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.
Computes and minimizes the power dissipation within given cellZones. In the absence of significant vi...
virtual void update_divDxDbMultiplier()
Update div(dx/db multiplier). Volume-based sensitivity term.
void update_dJdTMvar2()
Update field to be added to the second adjoint turbulence model PDE.
virtual void update_gradDxDbMultiplier()
Update grad(dx/db multiplier). Volume-based sensitivity term.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
Abstract base class for objective functions in incompressible flows.
virtual ~objectivePowerDissipation()=default
Destructor.
TypeName("powerDissipation")
Runtime type information.
Namespace for OpenFOAM.