scipy_minimizer.py 5.43 KB
Newer Older
Martin Reinecke's avatar
Martin Reinecke committed
1
2
3
4
5
6
7
8
9
10
11
12
13
# This program is free software: you can redistribute it and/or modify
# it under the terms of the GNU General Public License as published by
# the Free Software Foundation, either version 3 of the License, or
# (at your option) any later version.
#
# This program is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
# GNU General Public License for more details.
#
# You should have received a copy of the GNU General Public License
# along with this program.  If not, see <http://www.gnu.org/licenses/>.
#
Martin Reinecke's avatar
Martin Reinecke committed
14
# Copyright(C) 2013-2018 Max-Planck-Society
Martin Reinecke's avatar
Martin Reinecke committed
15
16
17
18
#
# NIFTy is being developed at the Max-Planck-Institut fuer Astrophysik
# and financially supported by the Studienstiftung des deutschen Volkes.

19
from __future__ import absolute_import, division, print_function
20
21

from .. import dobj
22
from ..compat import *
Martin Reinecke's avatar
Martin Reinecke committed
23
from ..field import Field
Martin Reinecke's avatar
Martin Reinecke committed
24
from ..logger import logger
25
from .iteration_controller import IterationController
26
from .minimizer import Minimizer
27
28


Martin Reinecke's avatar
Martin Reinecke committed
29
30
31
32
33
34
35
36
37
def _toNdarray(fld):
    return fld.to_global_data().reshape(-1)


def _toFlatNdarray(fld):
    return fld.val.flatten()


def _toField(arr, dom):
38
    return Field.from_global_data(dom, arr.reshape(dom.shape).copy())
Martin Reinecke's avatar
Martin Reinecke committed
39
40


41
42
43
44
45
46
class _MinHelper(object):
    def __init__(self, energy):
        self._energy = energy
        self._domain = energy.position.domain

    def _update(self, x):
Martin Reinecke's avatar
Martin Reinecke committed
47
48
        pos = _toField(x, self._domain)
        if (pos != self._energy.position).any():
49
            self._energy = self._energy.at(pos)
50
51
52
53
54
55
56

    def fun(self, x):
        self._update(x)
        return self._energy.value

    def jac(self, x):
        self._update(x)
Martin Reinecke's avatar
Martin Reinecke committed
57
        return _toFlatNdarray(self._energy.gradient)
58
59
60

    def hessp(self, x, p):
        self._update(x)
Martin Reinecke's avatar
Martin Reinecke committed
61
        res = self._energy.metric(_toField(p, self._domain))
Martin Reinecke's avatar
Martin Reinecke committed
62
        return _toFlatNdarray(res)
Martin Reinecke's avatar
Martin Reinecke committed
63
64
65
66
67
68
69
70
71
72
73
74
75


class ScipyMinimizer(Minimizer):
    """Scipy-based minimizer

    Parameters
    ----------
    method     : str
        The selected Scipy minimization method.
    options    : dictionary
        A set of custom options for the selected minimizer.
    """

76
    def __init__(self, method, options, need_hessp, bounds):
Martin Reinecke's avatar
Martin Reinecke committed
77
78
79
80
81
82
        super(ScipyMinimizer, self).__init__()
        if not dobj.is_numpy():
            raise NotImplementedError
        self._method = method
        self._options = options
        self._need_hessp = need_hessp
83
        self._bounds = bounds
Martin Reinecke's avatar
Martin Reinecke committed
84
85
86

    def __call__(self, energy):
        import scipy.optimize as opt
87
88
89
90
91
92
93
94
        hlp = _MinHelper(energy)
        energy = None  # drop handle, since we don't need it any more
        bounds = None
        if self._bounds is not None:
            if len(self._bounds) == 2:
                lo = self._bounds[0]
                hi = self._bounds[1]
                bounds = [(lo, hi)]*hlp._energy.position.size
Martin Reinecke's avatar
Martin Reinecke committed
95
            else:
96
97
98
99
100
101
                raise ValueError("unrecognized bounds")

        x = hlp._energy.position.val.flatten()
        hessp = hlp.hessp if self._need_hessp else None
        r = opt.minimize(hlp.fun, x, method=self._method, jac=hlp.jac,
                         hessp=hessp, options=self._options, bounds=bounds)
Martin Reinecke's avatar
Martin Reinecke committed
102
        if not r.success:
103
            logger.error("Problem in Scipy minimization: {}".format(r.message))
104
105
            return hlp._energy, IterationController.ERROR
        return hlp._energy, IterationController.CONVERGED
Martin Reinecke's avatar
Martin Reinecke committed
106
107


108
def NewtonCG(xtol, maxiter, disp=False):
Martin Reinecke's avatar
Martin Reinecke committed
109
110
111
112
113
114
    """Returns a ScipyMinimizer object carrying out the Newton-CG algorithm.

    See Also
    --------
    ScipyMinimizer
    """
Martin Reinecke's avatar
fix    
Martin Reinecke committed
115
    options = {"xtol": xtol, "maxiter": maxiter, "disp": disp}
116
    return ScipyMinimizer("Newton-CG", options, True, None)
Martin Reinecke's avatar
Martin Reinecke committed
117
118


119
def L_BFGS_B(ftol, gtol, maxiter, maxcor=10, disp=False, bounds=None):
Martin Reinecke's avatar
Martin Reinecke committed
120
121
122
123
124
125
    """Returns a ScipyMinimizer object carrying out the L-BFGS-B algorithm.

    See Also
    --------
    ScipyMinimizer
    """
126
    options = {"ftol": ftol, "gtol": gtol, "maxiter": maxiter,
Martin Reinecke's avatar
fix    
Martin Reinecke committed
127
               "maxcor": maxcor, "disp": disp}
128
    return ScipyMinimizer("L-BFGS-B", options, False, bounds)
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149


class ScipyCG(Minimizer):
    def __init__(self, tol, maxiter):
        super(ScipyCG, self).__init__()
        if not dobj.is_numpy():
            raise NotImplementedError
        self._tol = tol
        self._maxiter = maxiter

    def __call__(self, energy, preconditioner=None):
        from scipy.sparse.linalg import LinearOperator as scipy_linop, cg
        from .quadratic_energy import QuadraticEnergy
        if not isinstance(energy, QuadraticEnergy):
            raise ValueError("need a quadratic energy for CG")

        class mymatvec(object):
            def __init__(self, op):
                self._op = op

            def __call__(self, inp):
Martin Reinecke's avatar
Martin Reinecke committed
150
                return _toNdarray(self._op(_toField(inp, self._op.domain)))
151
152

        op = energy._A
Martin Reinecke's avatar
Martin Reinecke committed
153
154
        b = _toNdarray(energy._b)
        sx = _toNdarray(energy.position)
155
156
157
158
159
160
161
162
163
164
        sci_op = scipy_linop(shape=(op.domain.size, op.target.size),
                             matvec=mymatvec(op))
        prec_op = None
        if preconditioner is not None:
            prec_op = scipy_linop(shape=(op.domain.size, op.target.size),
                                  matvec=mymatvec(preconditioner))
        res, stat = cg(sci_op, b, x0=sx, tol=self._tol, M=prec_op,
                       maxiter=self._maxiter)
        stat = (IterationController.CONVERGED if stat >= 0 else
                IterationController.ERROR)
Martin Reinecke's avatar
Martin Reinecke committed
165
        return energy.at(_toField(res, op.domain)), stat