objectiveManager.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  Copyright (C) 2019 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
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 Class
29  Foam::objectiveManager
30 
31 Description
32  Class for managing objective functions.
33 
34 SourceFiles
35  objectiveManager.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef objectiveManager_H
40 #define objectiveManager_H
41 
42 #include "fvMesh.H"
43 #include "objective.H"
44 #include "runTimeSelectionTables.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 /*---------------------------------------------------------------------------*\
52  Class objectiveManager Declaration
53 \*---------------------------------------------------------------------------*/
54 
55 class objectiveManager
56 :
57  public regIOobject
58 {
59 protected:
60 
61  // Protected data
62 
63  const fvMesh& mesh_;
66  const word primalSolverName_;
69 
70 
71 private:
72 
73  // Private Member Functions
74 
75  //- No copy construct
76  objectiveManager(const objectiveManager&) = delete;
77 
78  //- No copy assignment
79  void operator=(const objectiveManager&) = delete;
80 
81 
82 public:
83 
84  TypeName("objectiveManager");
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 ~objectiveManager() = default;
100 
101 
102  // Member Functions
103 
104  virtual bool readDict(const dictionary& dict);
105 
106  //- Update objective function related values
108 
109  //- Update objective function related values
110  void update();
111 
112  //- Update contributions to adjoint if true, otherwise return nulls
113  void updateOrNullify();
114 
115  //- Increment integration times by the optimisation cycle time-span
116  void incrementIntegrationTimes(const scalar timeSpan);
117 
118  //- Print to screen
119  scalar print(bool negate = false);
120 
121  //- Should the objectives be written to file upon calling write()?
122  void setWrite(const bool shouldWrite);
123  //- Write objective function history
124  virtual bool writeObjectives();
125 
126  //- Write objective function history
127  virtual bool writeObjectives
128  (
129  const scalar weightedObjective,
130  const bool valid = true
131  );
132 
133  //- Call all functions required prior to the solution of the adjoint
134  //- equations
135  void updateAndWrite();
136 
137  //- Return reference to objective functions
139 
140  //- Return constant reference to objective functions
142 
143  //- Return name of adjointSolverManager
144  const word& adjointSolverManagerName() const;
145 
146  //- Return name of the adjointSolver
147  const word& adjointSolverName() const;
148 
149  //- Return name of the primalSolver
150  const word& primalSolverName() const;
151 
152  //- Check integration times for unsteady runs
153  void checkIntegrationTimes() const;
154 
155  //- Add contribution to adjoint momentum PDEs
156  virtual void addSource(fvVectorMatrix& matrix);
157 
158  //- Add contribution to a scalar adjoint PDEs
159  virtual void addSource(fvScalarMatrix& matrix);
160 
161 
162  // IO
163 
164  virtual bool writeData(Ostream&) const
165  {
166  return true;
167  }
168 };
169 
170 
171 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172 
173 } // End namespace Foam
174 
175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176 
177 #endif
178 
179 // ************************************************************************* //
void setWrite(const bool shouldWrite)
Should the objectives be written to file upon calling write()?
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
PtrList< objective > & getObjectiveFunctions()
Return reference to objective functions.
Class for managing objective functions.
void incrementIntegrationTimes(const scalar timeSpan)
Increment integration times by the optimisation cycle time-span.
const word & adjointSolverManagerName() const
Return name of adjointSolverManager.
void updateAndWrite()
Call all functions required prior to the solution of the adjoint equations.
PtrList< objective > objectives_
void updateNormalizationFactor()
Update objective function related values.
dynamicFvMesh & mesh
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual ~objectiveManager()=default
Destructor.
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
virtual bool writeObjectives()
Write objective function history.
virtual bool writeData(Ostream &) const
Pure virtual writeData function.
void negate(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1)
const word & primalSolverName() const
Return name of the primalSolver.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
scalar print(bool negate=false)
Print to screen.
void update()
Update objective function related values.
virtual void addSource(fvVectorMatrix &matrix)
Add contribution to adjoint momentum PDEs.
const word & adjointSolverName() const
Return name of the adjointSolver.
virtual bool readDict(const dictionary &dict)
void updateOrNullify()
Update contributions to adjoint if true, otherwise return nulls.
autoPtr< OFstream > weightedObjectiveFile_
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:66
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Macros to ease declaration of run-time selection tables.
TypeName("objectiveManager")
void checkIntegrationTimes() const
Check integration times for unsteady runs.
Namespace for OpenFOAM.