adjointSolver.C
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-2021 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
28\*---------------------------------------------------------------------------*/
29
30#include "adjointSolver.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
38}
39
40
41// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42
44(
45 fvMesh& mesh,
46 const word& managerType,
47 const dictionary& dict,
48 const word& primalSolverName
49)
50:
51 solver(mesh, managerType, dict),
52 primalSolverName_(primalSolverName),
53 objectiveManagerPtr_
54 (
56 (
57 mesh,
58 dict.subDict("objectives"),
59 solverName_,
60 primalSolverName
61 )
62 ),
63 sensitivities_(nullptr),
64 computeSensitivities_
65 (
66 dict.getOrDefault<bool>("computeSensitivities", true)
67 ),
68 isConstraint_(dict.getOrDefault<bool>("isConstraint", false))
69{
70 // Update objective-related quantities to get correct derivatives
71 // in case of continuation
72 objectiveManagerPtr_().update();
73}
74
75
76// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
77
79(
80 fvMesh& mesh,
81 const word& managerType,
82 const dictionary& dict,
83 const word& primalSolverName
84)
85{
86 const word solverType(dict.get<word>("type"));
87
88 auto* ctorPtr = adjointSolverConstructorTable(solverType);
89
90 if (!ctorPtr)
91 {
93 (
94 dict,
95 "adjointSolver",
96 solverType,
97 *adjointSolverConstructorTablePtr_
98 ) << exit(FatalIOError);
99 }
100
102 (
103 ctorPtr(mesh, managerType, dict, primalSolverName)
104 );
105}
106
107
108// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
109
111{
112 return mesh_.lookupObject<primalSolver>(primalSolverName_);
113}
114
115
117{
118 return
119 const_cast<primalSolver&>
120 (
121 mesh_.lookupObject<primalSolver>(primalSolverName_)
122 );
123}
124
125
127{
128 if (solver::readDict(dict))
129 {
130 computeSensitivities_ =
131 dict.getOrDefault<bool>("computeSensitivities", true);
132
133 objectiveManagerPtr_->readDict(dict.subDict("objectives"));
134
135 return true;
136 }
137
138 return false;
139}
140
141
143{
144 return objectiveManagerPtr_();
145}
146
147
149{
150 return objectiveManagerPtr_();
151}
152
153
155{
156 addProfiling(adjointSolver, "adjointSolver::postLoop");
157 computeObjectiveSensitivities();
158 // The solver dictionary has been already written after the termination
159 // of the adjoint loop. Force re-writing it to include the sensitivities
160 // as well
161 regIOobject::write(true);
162}
163
164
166{
167 return isConstraint_;
168}
169
170
172{
173 sensitivities_.clear();
174}
175
176
178{
179 // Does nothing in base
180}
181
182
184{
185 if (sensitivities_.valid())
186 {
187 sensitivities_().writeEntry("sensitivities", os);
188 }
189 return true;
190}
191
192
193// ************************************************************************* //
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Base class for adjoint solvers.
Definition: adjointSolver.H:60
virtual void clearSensitivities()
Clears the sensitivity field known by the adjoint solver.
autoPtr< objectiveManager > objectiveManagerPtr_
Object to manage objective functions.
Definition: adjointSolver.H:80
virtual bool writeData(Ostream &os) const
Write the sensitivity derivatives.
virtual bool readDict(const dictionary &dict)
const primalSolver & getPrimalSolver() const
virtual void updatePrimalBasedQuantities()
virtual bool isConstraint()
Is the solving referring to a constraint.
virtual void postLoop()
Functions to be called after loop.
const objectiveManager & getObjectiveManager() const
Return a const reference to the objective manager.
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
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
virtual bool write()
Write the output fields.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
class for managing incompressible objective functions.
Base class for primal solvers.
Definition: primalSolver.H:55
Base class for solution control classes.
Definition: solver.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
bool
Definition: EEqn.H:20
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
IOerror FatalIOError
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
#define addProfiling(name, descr)
Define profiling trigger with specified name and description string.
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict