adjointSimple.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::adjointSimple
30 
31 Description
32  Solution of the adjoint PDEs for incompressible, steady-state flows
33 
34  Reference:
35  \verbatim
36  For the development of the adjoint PDEs
37 
38  Papoutsis-Kiachagias, E. M., & Giannakoglou, K. C. (2014).
39  Continuous Adjoint Methods for Turbulent Flows, Applied to Shape
40  and Topology Optimization: Industrial Applications.
41  Archives of Computational Methods in Engineering, 23(2), 255-299.
42  http://doi.org/10.1007/s11831-014-9141-9
43  \endverbatim
44 
45 \*---------------------------------------------------------------------------*/
46 
47 #ifndef adjointSimple_H
48 #define adjointSimple_H
49 
51 #include "SIMPLEControl.H"
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 namespace Foam
56 {
57 
58 /*---------------------------------------------------------------------------*\
59  Class adjointSimple Declaration
60 \*---------------------------------------------------------------------------*/
61 
62 class adjointSimple
63 :
65 {
66 
67 private:
68 
69  // Private Member Functions
70 
71  //- No copy construct
72  adjointSimple(const adjointSimple&) = delete;
73 
74  //- No copy assignment
75  void operator=(const adjointSimple&) = delete;
76 
77 
78 protected:
79 
80  // Protected data
81 
82  //- Solver control
84 
85  //- Reference to incompressibleAdjointVars
86  // Used for convenience and to avoid repetitive dynamic_casts
87  // Same as getAdjointVars()
89 
90  //- Cumulative continuity error
91  scalar cumulativeContErr_;
92 
93 
94  // Protected Member Functions
95 
96  //- Allocate incompressibleAdjointVars and return reference to be used
97  //- for convenience in the rest of the class.
99 
100  //- Compute continuity errors
101  void continuityErrors();
102 
103  //- Accumulate the sensitivities integrand before calculating them
104  virtual void preCalculateSensitivities();
105 
106 
107 public:
108 
109  // Static Data Members
110 
111  //- Run-time type information
112  TypeName("adjointSimple");
113 
114 
115  // Constructors
116 
117  //- Construct from mesh and dictionary
119  (
120  fvMesh& mesh,
121  const word& managerType,
122  const dictionary& dict,
123  const word& primalSolverName,
124  const word& solverName
125  );
126 
127 
128  //- Destructor
129  virtual ~adjointSimple() = default;
130 
131 
132  // Member Functions
133 
134  // Evolution
135 
136  //- Execute one iteration of the solution algorithm
137  virtual void solveIter();
138 
139  //- Steps to be executed before each main SIMPLE iteration
140  virtual void preIter();
141 
142  //- The main SIMPLE iter
143  virtual void mainIter();
144 
145  //- Steps to be executed before each main SIMPLE iteration
146  virtual void postIter();
147 
148  //- Main control loop
149  virtual void solve();
150 
151  //- Looper (advances iters, time step)
152  virtual bool loop();
153 
154  //- Functions to be called before loop
155  virtual void preLoop();
156 
157 
158  // Source terms to be added to the adjoint equations
159 
160  //- Source terms for the momentum equations
161  virtual void addMomentumSource(fvVectorMatrix& matrix);
162 
163  //- Source terms for the continuity equation
164  virtual void addPressureSource(fvScalarMatrix& matrix);
165 
166 
167  //- Update primal based quantities
168  //- related to the objective functions.
169  // Also writes the objective function values to files.
170  // Written here and not in the postLoop function of the primal
171  // to make sure we don't pollute the objective files with
172  // objectives of non-converged linearSearch iterations
173  virtual void updatePrimalBasedQuantities();
174 
175  //- Add fvOptions for TopO
176  virtual void addTopOFvOptions() const;
177 
178  //- Write average iteration
179  virtual bool writeData(Ostream& os) const;
180 };
181 
182 
183 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
184 
185 } // End namespace Foam
186 
187 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 
189 #endif
190 
191 // ************************************************************************* //
virtual void mainIter()
The main SIMPLE iter.
virtual void addTopOFvOptions() const
Add fvOptions for TopO.
const word & primalSolverName() const
Return the primal solver name.
virtual void updatePrimalBasedQuantities()
Update primal based quantities related to the objective functions.
virtual ~adjointSimple()=default
Destructor.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual void preLoop()
Functions to be called before loop.
const fvMesh & mesh() const
Return the solver mesh.
Definition: solverI.H:24
autoPtr< SIMPLEControl > solverControl_
Solver control.
Definition: adjointSimple.H:82
TypeName("adjointSimple")
Run-time type information.
incompressibleAdjointVars & allocateVars()
Allocate incompressibleAdjointVars and return reference to be used for convenience in the rest of the...
Definition: adjointSimple.C:46
scalar cumulativeContErr_
Cumulative continuity error.
Definition: adjointSimple.H:95
Class including all adjoint fields for incompressible flows.
Base class for incompressibleAdjoint solvers.
virtual bool loop()
Looper (advances iters, time step)
incompressibleAdjointVars & adjointVars_
Reference to incompressibleAdjointVars.
Definition: adjointSimple.H:90
virtual void preIter()
Steps to be executed before each main SIMPLE iteration.
A class for handling words, derived from Foam::string.
Definition: word.H:63
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
const word & managerType() const
Return the manager type.
Definition: solverI.H:72
const dictionary & dict() const
Return the solver dictionary.
Definition: solverI.H:54
void continuityErrors()
Compute continuity errors.
Definition: adjointSimple.C:62
virtual void postIter()
Steps to be executed before each main SIMPLE iteration.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
OBJstream os(runTime.globalPath()/outputName)
virtual void solve()
Main control loop.
virtual void addPressureSource(fvScalarMatrix &matrix)
Source terms for the continuity equation.
const word & solverName() const
Return the solver name.
Definition: solverI.H:30
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Solution of the adjoint PDEs for incompressible, steady-state flows.
Definition: adjointSimple.H:55
virtual bool writeData(Ostream &os) const
Write average iteration.
virtual void preCalculateSensitivities()
Accumulate the sensitivities integrand before calculating them.
Definition: adjointSimple.C:81
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
virtual void solveIter()
Execute one iteration of the solution algorithm.
virtual void addMomentumSource(fvVectorMatrix &matrix)
Source terms for the momentum equations.
Namespace for OpenFOAM.