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-2020 PCOpt/NTUA
9  Copyright (C) 2013-2020 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 #include "incompressibleVars.H"
55 
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 
58 namespace Foam
59 {
60 
61 /*---------------------------------------------------------------------------*\
62  Class adjointSimple Declaration
63 \*---------------------------------------------------------------------------*/
64 
65 class adjointSimple
66 :
68 {
69 
70 private:
71 
72  // Private Member Functions
73 
74  //- No copy construct
75  adjointSimple(const adjointSimple&) = delete;
76 
77  //- No copy assignment
78  void operator=(const adjointSimple&) = delete;
79 
80 
81 protected:
82 
83  // Protected data
84 
85  //- Solver control
87 
88  //- Reference to incompressibleAdjointVars
89  // Used for convenience and to avoid repetitive dynamic_casts
90  // Same as getAdjointVars()
92 
93  //- Cumulative continuity error
94  scalar cumulativeContErr_;
95 
96  //- Sensitivity Derivatives engine
98 
99 
100  // Protected Member Functions
101 
102  //- Allocate incompressibleAdjointVars and return reference to be used
103  //- for convenience in the rest of the class.
105 
106  //- In case variable names are different than the base ones,
107  //- add extra schemes and relaxation factors to the appropriate dicts
108  // Note: 160813: Changes in the baseline solution and fvSchemes classes
109  // have to be made in order to add schemes in the dict at run time.
110  // Not supported for now
111  void addExtraSchemes();
112 
113  //- Compute continuity errors
114  void continuityErrors();
115 
116 
117 public:
118 
119  // Static Data Members
120 
121  //- Run-time type information
122  TypeName("adjointSimple");
123 
124 
125  // Constructors
126 
127  //- Construct from mesh and dictionary
129  (
130  fvMesh& mesh,
131  const word& managerType,
132  const dictionary& dict,
133  const word& primalSolverName
134  );
135 
136 
137  //- Destructor
138  virtual ~adjointSimple() = default;
139 
140 
141  // Member Functions
142 
143  //- Read dict if updated
144  virtual bool readDict(const dictionary& dict);
145 
146 
147  // Evolution
148 
149  //- Execute one iteration of the solution algorithm
150  virtual void solveIter();
151 
152  //- Steps to be executed before each main SIMPLE iteration
153  virtual void preIter();
154 
155  //- The main SIMPLE iter
156  virtual void mainIter();
157 
158  //- Steps to be executed before each main SIMPLE iteration
159  virtual void postIter();
160 
161  //- Main control loop
162  virtual void solve();
163 
164  //- Looper (advances iters, time step)
165  virtual bool loop();
166 
167  //- Functions to be called before loop
168  virtual void preLoop();
169 
170  //- Compute sensitivities of the underlaying objectives
171  virtual void computeObjectiveSensitivities();
172 
173  //- Grab a reference to the computed sensitivities
174  virtual const scalarField& getObjectiveSensitivities();
175 
176  //- Clears the sensitivity field known by the adjoint solver
177  //- and zeros sensitivities constituents
178  virtual void clearSensitivities();
179 
180  //- Return the base sensitivity object
181  virtual sensitivity& getSensitivityBase();
182 
183 
184  // Source terms to be added to the adjoint equations
185 
186  //- Source terms for the momentum equations
187  virtual void addMomentumSource(fvVectorMatrix& matrix);
188 
189  //- Source terms for the continuity equation
190  virtual void addPressureSource(fvScalarMatrix& matrix);
191 
192  //- Update primal based quantities
193  //- related to the objective functions.
194  // Also writes the objective function values to files.
195  // Written here and not in the postLoop function of the primal
196  // to make sure we don't pollute the objective files with
197  // objectives of non-converged linearSearch iterations
198  virtual void updatePrimalBasedQuantities();
199 
200  //- Write average iteration
201  virtual bool writeData(Ostream& os) const;
202 };
203 
204 
205 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
206 
207 } // End namespace Foam
208 
209 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
210 
211 #endif
212 
213 // ************************************************************************* //
Foam::adjointSimple::continuityErrors
void continuityErrors()
Compute continuity errors.
Definition: adjointSimple.C:83
Foam::adjointSimple::clearSensitivities
virtual void clearSensitivities()
Definition: adjointSimple.C:386
Foam::incompressibleAdjointSolver
Base class for incompressibleAdjoint solvers.
Definition: incompressibleAdjointSolver.H:53
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::adjointSimple::solveIter
virtual void solveIter()
Execute one iteration of the solution algorithm.
Definition: adjointSimple.C:182
Foam::adjointSolver::primalSolverName
const word & primalSolverName() const
Return the primal solver name.
Definition: adjointSolver.H:152
incompressibleVars.H
Foam::adjointSimple::computeObjectiveSensitivities
virtual void computeObjectiveSensitivities()
Compute sensitivities of the underlaying objectives.
Definition: adjointSimple.C:356
Foam::adjointSimple::addMomentumSource
virtual void addMomentumSource(fvVectorMatrix &matrix)
Source terms for the momentum equations.
Definition: adjointSimple.C:412
Foam::adjointSimple::adjointVars_
incompressibleAdjointVars & adjointVars_
Reference to incompressibleAdjointVars.
Definition: adjointSimple.H:90
Foam::adjointSimple::preLoop
virtual void preLoop()
Functions to be called before loop.
Definition: adjointSimple.C:349
Foam::incompressibleAdjointVars
Class including all adjoint fields for incompressible flows.
Definition: incompressibleAdjointVars.H:52
Foam::adjointSimple::preIter
virtual void preIter()
Steps to be executed before each main SIMPLE iteration.
Definition: adjointSimple.C:190
incompressibleAdjointSolver.H
Foam::adjointSimple::adjointSensitivity_
autoPtr< incompressible::adjointSensitivity > adjointSensitivity_
Sensitivity Derivatives engine.
Definition: adjointSimple.H:96
Foam::Field< scalar >
Foam::adjointSimple::addPressureSource
virtual void addPressureSource(fvScalarMatrix &matrix)
Source terms for the continuity equation.
Definition: adjointSimple.C:418
Foam::adjointSimple::getSensitivityBase
virtual sensitivity & getSensitivityBase()
Return the base sensitivity object.
Definition: adjointSimple.C:396
SIMPLEControl.H
Foam::adjointSimple::~adjointSimple
virtual ~adjointSimple()=default
Destructor.
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::adjointSimple::solve
virtual void solve()
Main control loop.
Definition: adjointSimple.C:330
Foam::adjointSimple
Solution of the adjoint PDEs for incompressible, steady-state flows.
Definition: adjointSimple.H:64
Foam::adjointSimple::mainIter
virtual void mainIter()
The main SIMPLE iter.
Definition: adjointSimple.C:196
incompressibleAdjointVars.H
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::solver::mesh
const fvMesh & mesh() const
Return the solver mesh.
Definition: solver.C:96
Foam::adjointSimple::readDict
virtual bool readDict(const dictionary &dict)
Read dict if updated.
Definition: adjointSimple.C:160
adjointSensitivityIncompressible.H
Foam::adjointSimple::updatePrimalBasedQuantities
virtual void updatePrimalBasedQuantities()
Definition: adjointSimple.C:424
Foam::adjointSimple::postIter
virtual void postIter()
Steps to be executed before each main SIMPLE iteration.
Definition: adjointSimple.C:318
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::adjointSimple::writeData
virtual bool writeData(Ostream &os) const
Write average iteration.
Definition: adjointSimple.C:433
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::adjointSimple::TypeName
TypeName("adjointSimple")
Run-time type information.
Foam::adjointSimple::getObjectiveSensitivities
virtual const scalarField & getObjectiveSensitivities()
Grab a reference to the computed sensitivities.
Definition: adjointSimple.C:375
Foam::sensitivity
Abstract base class for adjoint sensitivities.
Definition: sensitivity.H:63
Foam::adjointSimple::addExtraSchemes
void addExtraSchemes()
Definition: adjointSimple.C:69
Foam::solver::dict
virtual const dictionary & dict() const
Return the solver dictionary.
Definition: solver.C:113
Foam::adjointSimple::loop
virtual bool loop()
Looper (advances iters, time step)
Definition: adjointSimple.C:343
Foam::adjointSimple::solverControl_
autoPtr< SIMPLEControl > solverControl_
Solver control.
Definition: adjointSimple.H:85
Foam::adjointSimple::cumulativeContErr_
scalar cumulativeContErr_
Cumulative continuity error.
Definition: adjointSimple.H:93
Foam::adjointSimple::allocateVars
incompressibleAdjointVars & allocateVars()
Definition: adjointSimple.C:53