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-2021 PCOpt/NTUA
9 Copyright (C) 2013-2021 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
28Class
29 Foam::objectiveManager
30
31Description
32 class for managing incompressible objective functions.
33
34SourceFiles
35 objectiveManager.C
36
37\*---------------------------------------------------------------------------*/
38
39#ifndef objectiveManager_H
40#define objectiveManager_H
41
42#include "fvMesh.H"
43#include "objective.H"
45
46// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47
48namespace Foam
49{
50
51/*---------------------------------------------------------------------------*\
52 Class objectiveManager Declaration
53\*---------------------------------------------------------------------------*/
56:
57 public regIOobject
58{
59protected:
60
61 // Protected data
63 const fvMesh& mesh_;
69
70
71private:
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
82public:
84 TypeName("objectiveManager");
85
86 // Declare run-time constructor selection table
89 (
90 autoPtr,
93 (
94 const fvMesh& mesh,
95 const dictionary& dict,
98 ),
100 );
101
102 // Constructors
103
104 //- Construct from components
106 (
107 const fvMesh& mesh,
108 const dictionary& dict,
109 const word& adjointSolverName,
111 );
112
113 // Selectors
114
115 //- Return a reference to the selected turbulence model
117 (
118 const fvMesh& mesh,
119 const dictionary& dict,
120 const word& adjointSolverName,
122 );
123
124
125 //- Destructor
126 virtual ~objectiveManager() = default;
127
128
129 // Member Functions
130
131 virtual bool readDict(const dictionary& dict);
132
133 //- Update objective function related values
135
136 //- Update objective function related values
137 void update();
138
139 //- Update contributions to adjoint if true, otherwise return nulls
140 void updateOrNullify();
141
142 //- Increment integration times by the optimisation cycle time-span
143 void incrementIntegrationTimes(const scalar timeSpan);
144
145 //- Print to screen
146 scalar print();
147
148 //- Write objective function history
149 virtual bool writeObjectives
150 (
151 const scalar weightedObjective,
152 const bool valid = true
153 );
154
155 //- Call all functions required prior to the solution of the adjoint
156 //- equations
157 void updateAndWrite();
158
159 //- Return reference to objective functions
161
162 //- Return constant reference to objective functions
164
165 //- Return name of adjointSolverManager
166 const word& adjointSolverManagerName() const;
167
168 //- Return name of the adjointSolver
169 const word& adjointSolverName() const;
170
171 //- Return name of the primalSolver
172 const word& primalSolverName() const;
173
174 //- Check integration times for unsteady runs
175 void checkIntegrationTimes() const;
176
177 //- Add contribution to adjoint momentum PDEs
178 virtual void addUaEqnSource(fvVectorMatrix& UaEqn) = 0;
179
180 //- Add contribution to adjoint momentum PDEs
181 virtual void addPaEqnSource(fvScalarMatrix& paEqn) = 0;
182
183 //- Add contribution to first adjoint turbulence model PDE
184 virtual void addTMEqn1Source(fvScalarMatrix& adjTMEqn1) = 0;
185
186 //- Add contribution to second adjoint turbulence model PDE
187 virtual void addTMEqn2Source(fvScalarMatrix& adjTMEqn2) = 0;
188
189
190 // IO
192 virtual bool writeData(Ostream&) const
193 {
194 return true;
195 }
196};
197
198
199// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
200
201} // End namespace Foam
202
203// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204
205#endif
206
207// ************************************************************************* //
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
class for managing incompressible objective functions.
virtual void addTMEqn2Source(fvScalarMatrix &adjTMEqn2)=0
Add contribution to second adjoint turbulence model PDE.
autoPtr< OFstream > weigthedObjectiveFile_
virtual bool writeObjectives(const scalar weightedObjective, const bool valid=true)
Write objective function history.
PtrList< objective > & getObjectiveFunctions()
Return reference to objective functions.
const word & adjointSolverName() const
Return name of the adjointSolver.
scalar print()
Print to screen.
declareRunTimeSelectionTable(autoPtr, objectiveManager, dictionary,(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName),(mesh, dict, adjointSolverName, primalSolverName))
virtual void addUaEqnSource(fvVectorMatrix &UaEqn)=0
Add contribution to adjoint momentum PDEs.
virtual void addPaEqnSource(fvScalarMatrix &paEqn)=0
Add contribution to adjoint momentum PDEs.
static autoPtr< objectiveManager > New(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName)
Return a reference to the selected turbulence model.
virtual ~objectiveManager()=default
Destructor.
virtual bool readDict(const dictionary &dict)
const dictionary & dict_
void checkIntegrationTimes() const
Check integration times for unsteady runs.
PtrList< objective > objectives_
TypeName("objectiveManager")
virtual bool writeData(Ostream &) const
Pure virtual writeData function.
const word & adjointSolverManagerName() const
Return name of adjointSolverManager.
virtual void addTMEqn1Source(fvScalarMatrix &adjTMEqn1)=0
Add contribution to first adjoint turbulence model PDE.
void updateOrNullify()
Update contributions to adjoint if true, otherwise return nulls.
void incrementIntegrationTimes(const scalar timeSpan)
Increment integration times by the optimisation cycle time-span.
void update()
Update objective function related values.
void updateNormalizationFactor()
Update objective function related values.
const word & primalSolverName() const
Return name of the primalSolver.
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:76
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
Namespace for OpenFOAM.
fvScalarMatrix paEqn(fvm::d2dt2(pa) - sqr(c0) *fvc::laplacian(pa))
Macros to ease declaration of run-time selection tables.
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes)
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73