1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218
|
// Copyright ©2014 The Gonum Authors. All rights reserved.
// Use of this source code is governed by a BSD-style
// license that can be found in the LICENSE file.
package optimize
import (
"math"
"gonum.org/v1/gonum/floats"
)
// LinesearchMethod represents an abstract optimization method in which a
// function is optimized through successive line search optimizations.
type LinesearchMethod struct {
// NextDirectioner specifies the search direction of each linesearch.
NextDirectioner NextDirectioner
// Linesearcher performs a linesearch along the search direction.
Linesearcher Linesearcher
x []float64 // Starting point for the current iteration.
dir []float64 // Search direction for the current iteration.
first bool // Indicator of the first iteration.
nextMajor bool // Indicates that MajorIteration must be commanded at the next call to Iterate.
eval Operation // Indicator of valid fields in Location.
lastStep float64 // Step taken from x in the previous call to Iterate.
lastOp Operation // Operation returned from the previous call to Iterate.
}
func (ls *LinesearchMethod) Init(loc *Location) (Operation, error) {
if loc.Gradient == nil {
panic("linesearch: gradient is nil")
}
dim := len(loc.X)
ls.x = resize(ls.x, dim)
ls.dir = resize(ls.dir, dim)
ls.first = true
ls.nextMajor = false
// Indicate that all fields of loc are valid.
ls.eval = FuncEvaluation | GradEvaluation
if loc.Hessian != nil {
ls.eval |= HessEvaluation
}
ls.lastStep = math.NaN()
ls.lastOp = NoOperation
return ls.initNextLinesearch(loc)
}
func (ls *LinesearchMethod) Iterate(loc *Location) (Operation, error) {
switch ls.lastOp {
case NoOperation:
// TODO(vladimir-ch): Either Init has not been called, or the caller is
// trying to resume the optimization run after Iterate previously
// returned with an error. Decide what is the proper thing to do. See also #125.
case MajorIteration:
// The previous updated location did not converge the full
// optimization. Initialize a new Linesearch.
return ls.initNextLinesearch(loc)
default:
// Update the indicator of valid fields of loc.
ls.eval |= ls.lastOp
if ls.nextMajor {
ls.nextMajor = false
// Linesearcher previously finished, and the invalid fields of loc
// have now been validated. Announce MajorIteration.
ls.lastOp = MajorIteration
return ls.lastOp, nil
}
}
// Continue the linesearch.
f := math.NaN()
if ls.eval&FuncEvaluation != 0 {
f = loc.F
}
projGrad := math.NaN()
if ls.eval&GradEvaluation != 0 {
projGrad = floats.Dot(loc.Gradient, ls.dir)
}
op, step, err := ls.Linesearcher.Iterate(f, projGrad)
if err != nil {
return ls.error(err)
}
switch op {
case MajorIteration:
// Linesearch has been finished.
ls.lastOp = complementEval(loc, ls.eval)
if ls.lastOp == NoOperation {
// loc is complete, MajorIteration can be declared directly.
ls.lastOp = MajorIteration
} else {
// Declare MajorIteration on the next call to Iterate.
ls.nextMajor = true
}
case FuncEvaluation, GradEvaluation, FuncEvaluation | GradEvaluation:
if step != ls.lastStep {
// We are moving to a new location, and not, say, evaluating extra
// information at the current location.
// Compute the next evaluation point and store it in loc.X.
floats.AddScaledTo(loc.X, ls.x, step, ls.dir)
if floats.Equal(ls.x, loc.X) {
// Step size has become so small that the next evaluation point is
// indistinguishable from the starting point for the current
// iteration due to rounding errors.
return ls.error(ErrNoProgress)
}
ls.lastStep = step
ls.eval = NoOperation // Indicate all invalid fields of loc.
}
ls.lastOp = op
default:
panic("linesearch: Linesearcher returned invalid operation")
}
return ls.lastOp, nil
}
func (ls *LinesearchMethod) error(err error) (Operation, error) {
ls.lastOp = NoOperation
return ls.lastOp, err
}
// initNextLinesearch initializes the next linesearch using the previous
// complete location stored in loc. It fills loc.X and returns an evaluation
// to be performed at loc.X.
func (ls *LinesearchMethod) initNextLinesearch(loc *Location) (Operation, error) {
copy(ls.x, loc.X)
var step float64
if ls.first {
ls.first = false
step = ls.NextDirectioner.InitDirection(loc, ls.dir)
} else {
step = ls.NextDirectioner.NextDirection(loc, ls.dir)
}
projGrad := floats.Dot(loc.Gradient, ls.dir)
if projGrad >= 0 {
return ls.error(ErrNonDescentDirection)
}
op := ls.Linesearcher.Init(loc.F, projGrad, step)
switch op {
case FuncEvaluation, GradEvaluation, FuncEvaluation | GradEvaluation:
default:
panic("linesearch: Linesearcher returned invalid operation")
}
floats.AddScaledTo(loc.X, ls.x, step, ls.dir)
if floats.Equal(ls.x, loc.X) {
// Step size is so small that the next evaluation point is
// indistinguishable from the starting point for the current iteration
// due to rounding errors.
return ls.error(ErrNoProgress)
}
ls.lastStep = step
ls.eval = NoOperation // Invalidate all fields of loc.
ls.lastOp = op
return ls.lastOp, nil
}
// ArmijoConditionMet returns true if the Armijo condition (aka sufficient
// decrease) has been met. Under normal conditions, the following should be
// true, though this is not enforced:
// - initGrad < 0
// - step > 0
// - 0 < decrease < 1
func ArmijoConditionMet(currObj, initObj, initGrad, step, decrease float64) bool {
return currObj <= initObj+decrease*step*initGrad
}
// StrongWolfeConditionsMet returns true if the strong Wolfe conditions have been met.
// The strong Wolfe conditions ensure sufficient decrease in the function
// value, and sufficient decrease in the magnitude of the projected gradient.
// Under normal conditions, the following should be true, though this is not
// enforced:
// - initGrad < 0
// - step > 0
// - 0 <= decrease < curvature < 1
func StrongWolfeConditionsMet(currObj, currGrad, initObj, initGrad, step, decrease, curvature float64) bool {
if currObj > initObj+decrease*step*initGrad {
return false
}
return math.Abs(currGrad) < curvature*math.Abs(initGrad)
}
// WeakWolfeConditionsMet returns true if the weak Wolfe conditions have been met.
// The weak Wolfe conditions ensure sufficient decrease in the function value,
// and sufficient decrease in the value of the projected gradient. Under normal
// conditions, the following should be true, though this is not enforced:
// - initGrad < 0
// - step > 0
// - 0 <= decrease < curvature< 1
func WeakWolfeConditionsMet(currObj, currGrad, initObj, initGrad, step, decrease, curvature float64) bool {
if currObj > initObj+decrease*step*initGrad {
return false
}
return currGrad >= curvature*initGrad
}
|