SR1.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 -------------------------------------------------------------------------------
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 \*---------------------------------------------------------------------------*/
29 
30 #include "SR1.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(SR1, 0);
39  (
40  updateMethod,
41  SR1,
43  );
44 }
45 
46 
47 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
48 
50 {
51  // Vectors needed to construct the inverse HessianInv matrix
56 
57  scalarField temp(s - rightMult(Hessian_(), y));
58 
59  // Construct the inverse HessianInv
60  scalar tempMag = sqrt(globalSum(sqr(temp)));
61  scalar yMag = sqrt(globalSum(sqr(y)));
62  scalar HessYMag = sqrt(globalSum(sqr(rightMult(Hessian_(), y))));
63 
64  // Stability check
65  if (tempMag > ratioThreshold_ * yMag * HessYMag)
66  {
67  Hessian_() += (scalar(1)/(globalSum(temp*y)))*outerProd(temp, temp);
68  }
69  else
70  {
72  << "Denominator of update too small. Keeping old Hessian" << endl;
73  }
74 }
75 
76 
77 void Foam::SR1::update()
78 {
79  // In the first few iterations, use steepest descent but update the Hessian
80  // matrix
81  if (counter_ < nSteepestDescent_)
82  {
83  Info<< "Using steepest descent to update design variables" << endl;
84  for (const label varI : activeDesignVars_)
85  {
86  correction_[varI] = -eta_*objectiveDerivatives_[varI];
87  }
88  }
89  else
90  {
91  scalarField activeDerivs(activeDesignVars_.size(), Zero);
92  activeDerivs.map(objectiveDerivatives_, activeDesignVars_);
93  scalarField activeCorrection
94  (
95  -etaHessian_*rightMult(Hessian_(), activeDerivs)
96  );
97 
98  // Transfer correction to the global list
99  correction_ = Zero;
100  forAll(activeDesignVars_, varI)
101  {
102  correction_[activeDesignVars_[varI]] = activeCorrection[varI];
103  }
104  }
105 
106  // Store fields for the next iteration
107  derivativesOld_ = objectiveDerivatives_;
108  correctionOld_ = correction_;
109 }
110 
111 
112 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
113 
114 Foam::SR1::SR1
115 (
116  const fvMesh& mesh,
117  const dictionary& dict,
118  autoPtr<designVariables>& designVars,
119  const label nConstraints,
120  const word& type
121 )
122 :
123  quasiNewton(mesh, dict, designVars, nConstraints, type),
124  ratioThreshold_
125  (
126  coeffsDict(type).getOrDefault<scalar>("ratioThreshold", 1e-08)
127  )
128 {
129  allocateHessian();
130 }
131 
132 
133 // ************************************************************************* //
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const scalarField rightMult(const SquareMatrix< scalar > &, const scalarField &)
Definition: updateMethod.C:65
scalar ratioThreshold_
For stability check.
Definition: SR1.H:61
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
scalar y
const labelList & activeDesignVars_
Map to active design variables.
Definition: updateMethod.H:75
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
Base class for quasi-Newton methods.
Definition: quasiNewton.H:49
A class for handling words, derived from Foam::string.
Definition: word.H:63
autoPtr< SquareMatrix< scalar > > Hessian_
The Hessian or its inverse, depending on the deriving class.
Definition: quasiNewton.H:78
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
scalarField derivativesOld_
The previous derivatives.
Definition: quasiNewton.H:83
void map(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 map from the given field
Definition: Field.C:303
virtual void update()
Update design variables.
Definition: SR1.C:70
scalar globalSum(const scalarField &field)
Compute either global or local sum, based on globalSum flag.
Definition: updateMethod.C:168
void allocateHessian()
Allocate the Hessian matrix.
Definition: quasiNewton.C:35
defineTypeNameAndDebug(combustionModel, 0)
scalarField objectiveDerivatives_
Derivatives of the objective functions.
Definition: updateMethod.H:80
#define WarningInFunction
Report a warning using Foam::Warning.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
SquareMatrix< scalar > outerProd(const scalarField &, const scalarField &)
Definition: updateMethod.C:92
virtual void updateHessian()
Update approximation of the inverse Hessian.
Definition: SR1.C:42
messageStream Info
Information stream (stdout output on master, null elsewhere)
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
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.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
scalarField correctionOld_
The previous correction.
Definition: quasiNewton.H:88