Loading...
Searching...
No Matches
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-2021 PCOpt/NTUA
9 Copyright (C) 2013-2021 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
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
39 (
43 );
44}
45
46
47// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48
49Foam::conjugateGradient::conjugateGradient
50(
51 const fvMesh& mesh,
52 const dictionary& dict,
53 autoPtr<designVariables>& designVars,
54 const label nConstraints,
55 const word& type
56)
57:
58 updateMethod(mesh, dict, designVars, nConstraints, type),
59
60 dxOld_(readOrZeroField("dxOld", activeDesignVars_.size())),
61 sOld_(readOrZeroField("sOld", activeDesignVars_.size())),
62 betaType_(coeffsDict(type).getOrDefault<word>("betaType", "FletcherReeves"))
63{
64 // Check if beta type is valid
65 if
66 (
67 !(betaType_ == "FletcherReeves")
68 && !(betaType_ == "PolakRibiere")
69 && !(betaType_ == "PolakRibiereRestarted")
70 )
71 {
73 << "Invalid betaType " << betaType_ << ". Valid options are "
74 << "FletcherReeves, PolakRibiere, PolakRibiereRestarted"
75 << nl << nl
77 }
78}
79
80
81// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
82
84{
85 if (counter_ == 0)
86 {
87 Info<< "Using steepest descent for the first iteration" << endl;
88 for (const label varI : activeDesignVars_)
89 {
90 correction_[varI] = -eta_*objectiveDerivatives_[varI];
91 }
92
93 dxOld_.map(-objectiveDerivatives_, activeDesignVars_);
94 sOld_ = dxOld_;
95 }
96 else
97 {
98 scalarField dx(activeDesignVars_.size(), Zero);
99 dx.map(-objectiveDerivatives_, activeDesignVars_);
100
101 scalar beta(Zero);
102 if (betaType_ == "FletcherReeves")
103 {
104 beta = globalSum(dx*dx)/globalSum(dxOld_ * dxOld_);
105 }
106 else if (betaType_ == "PolakRibiere")
107 {
108 beta = globalSum(dx*(dx - dxOld_))/globalSum(dxOld_ * dxOld_);
109 }
110 else if (betaType_ == "PolakRibiereRestarted")
111 {
112 beta =
113 max
114 (
115 scalar(0),
116 globalSum(dx*(dx - dxOld_))/globalSum(dxOld_ * dxOld_)
117 );
118 if (beta == scalar(0))
119 {
120 Info<< "Computed negative beta. Resetting to zero" << endl;
121 }
122 }
123
124 scalarField s(dx + beta*sOld_);
125
126 correction_ = Zero;
127 forAll(activeDesignVars_, varI)
128 {
129 correction_[activeDesignVars_[varI]] = eta_*s[varI];
130 }
131
132 // Store fields for the next iteration
133 dxOld_ = dx;
135 }
136
137 ++counter_;
138}
139
140
142(
143 const scalarField& oldCorrection
144)
146 sOld_.map(oldCorrection, activeDesignVars_);
147 sOld_ /= eta_;
148 correction_ = oldCorrection;
149}
150
151
153{
154 dxOld_.writeEntry("dxOld", os);
155 sOld_.writeEntry("sOld", os);
156
158}
159
160
161// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
label size() const noexcept
The number of elements in list.
Definition DLListBase.H:194
void map(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 map from the given field
Definition Field.C:302
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
The Conjugate Gradient formula.
void computeCorrection()
Compute design variables correction.
virtual bool writeData(Ostream &os) const
Write useful quantities to files.
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:133
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Abstract base class for optimisation methods.
scalarField objectiveDerivatives_
Derivatives of the objective functions.
scalarField correction_
Design variables correction.
dictionary coeffsDict(const word &type) const
Return optional dictionary with parameters specific to each method.
scalar globalSum(const scalarField &field)
Compute either global or local sum, based on globalSum flag.
virtual bool writeData(Ostream &os) const
Write continuation data under uniform.
tmp< scalarField > readOrZeroField(const word &name, const label size)
Helper function to either read a scalarField of certain size from a dictionary, or construct a zero f...
const labelList & activeDesignVars_
Map to active design variables.
label counter_
Optimisation cycle count.
label nConstraints() const
Get the number of constraints.
scalar eta_
Step multiplying the correction.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
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:40
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:519
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299