conjugateGradient.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-2019 PCOpt/NTUA
9 Copyright (C) 2013-2019 FOSS GP
10 Copyright (C) 2019-2020 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 "conjugateGradient.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
39 (
43 );
44}
45
46
47// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48
50{
51 // Set active design variables, if necessary
53 {
55 }
56
57 // Allocate old fields
60}
61
62
64{
65 if (optMethodIODict_.headerOk())
66 {
67 optMethodIODict_.readEntry("dxOld", dxOld_);
68 optMethodIODict_.readEntry("sOld", sOld_);
69 optMethodIODict_.readEntry("counter", counter_);
70 optMethodIODict_.readEntry("eta", eta_);
71
72 label nDVs = optMethodIODict_.get<label>("nDVs");
73 correction_ = scalarField(nDVs, Zero);
74
75 if (activeDesignVars_.empty())
76 {
77 activeDesignVars_ = identity(nDVs);
78 }
79 }
80}
81
82
83// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
84
86(
87 const fvMesh& mesh,
88 const dictionary& dict
89)
90:
92
93 activeDesignVars_(0),
94 dxOld_(0),
95 sOld_(0),
96 counter_(0),
97 betaType_
98 (
99 coeffsDict().getOrDefault<word>("betaType", "FletcherReeves")
100 )
101{
102 if
103 (
104 !coeffsDict().readIfPresent("activeDesignVariables", activeDesignVars_)
105 )
106 {
107 // If not, all available design variables will be used.
108 // Number is not know at the moment
109 Info<< "\t Did not find explicit definition of active design variables. "
110 << "Treating all available ones as active " << endl;
111 }
112
113 // Check if beta type is valid
114 if
115 (
116 !(betaType_ == "FletcherReeves")
117 && !(betaType_ == "PolakRibiere")
118 && !(betaType_ == "PolakRibiereRestarted")
119 )
120 {
122 << "Invalid betaType " << betaType_ << ". Valid options are "
123 << "FletcherReeves, PolakRibiere, PolakRibiereRestarted"
124 << nl << nl
125 << exit(FatalError);
126 }
127
128 // Read old dx and s, if present
129 readFromDict();
130}
131
132
133// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
134
136{
137 if (counter_ == 0)
138 {
139 allocateFields();
140
141 Info<< "Using steepest descent for the first iteration" << endl;
142 correction_ = -eta_*objectiveDerivatives_;
143
144 dxOld_.map(-objectiveDerivatives_, activeDesignVars_);
145 sOld_ = dxOld_;
146 }
147 else
148 {
149 scalarField dx = scalarField(activeDesignVars_.size(), Zero);
150 dx.map(-objectiveDerivatives_, activeDesignVars_);
151
152 scalar beta(Zero);
153 if (betaType_ == "FletcherReeves")
154 {
155 beta = globalSum(dx*dx)/globalSum(dxOld_ * dxOld_);
156 }
157 else if (betaType_ == "PolakRibiere")
158 {
159 beta = globalSum(dx*(dx - dxOld_))/globalSum(dxOld_ * dxOld_);
160 }
161 else if (betaType_ == "PolakRibiereRestarted")
162 {
163 beta =
164 max
165 (
166 scalar(0),
167 globalSum(dx*(dx - dxOld_))/globalSum(dxOld_ * dxOld_)
168 );
169 if (beta == scalar(0))
170 {
171 Info<< "Computed negative beta. Resetting to zero" << endl;
172 }
173 }
174
175 scalarField s(dx + beta*sOld_);
176
177 correction_ = Zero;
178 forAll(activeDesignVars_, varI)
179 {
180 correction_[activeDesignVars_[varI]] = eta_*s[varI];
181 }
182
183 // Store fields for the next iteration
184 dxOld_ = dx;
185 sOld_ = s;
186 }
187
188 ++counter_;
189}
190
191
193(
194 const scalarField& oldCorrection
195)
196{
197 sOld_.map(oldCorrection, activeDesignVars_);
198 sOld_ /= eta_;
199 correction_ = oldCorrection;
200}
201
202
204{
205 optMethodIODict_.add<scalarField>("dxOld", dxOld_, true);
206 optMethodIODict_.add<scalarField>("sOld", sOld_, true);
207 optMethodIODict_.add<label>("counter", counter_, true);
208 optMethodIODict_.add<label>("nDVs", objectiveDerivatives_.size(), true);
209
211}
212
213
214// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
void map(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 map from the given field
Definition: Field.C:240
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
The Conjugate Gradient formula.
void computeCorrection()
Compute design variables correction.
void readFromDict()
Read old info from dict.
void allocateFields()
Allocate matrices in the first optimisation cycle.
virtual void write()
Write old info to dict.
virtual void updateOldCorrection(const scalarField &oldCorrection)
Update old correction. For use when eta has been changed externally.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Abstract base class for optimisation methods.
Definition: updateMethod.H:55
scalarField objectiveDerivatives_
Derivatives of the objective functions.
Definition: updateMethod.H:68
dictionary coeffsDict()
Return optional dictionary with parameters specific to each method.
Definition: updateMethod.C:254
virtual void write()
Write useful quantities to files.
Definition: updateMethod.C:398
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
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:38
messageStream Info
Information stream (stdout output on master, null elsewhere)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
propsDict readIfPresent("fields", acceptFields)