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 20
from __future__ import absolute_import, division, print_function
from ..compat import *
Martin Reinecke's avatar
Martin Reinecke committed
21 22 23
from .minimizer import Minimizer
from ..field import Field
from .. import dobj
Martin Reinecke's avatar
Martin Reinecke committed
24
from ..logger import logger
25 26 27
from .iteration_controller import IterationController


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


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


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


40 41 42 43 44 45
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
46 47
        pos = _toField(x, self._domain)
        if (pos != self._energy.position).any():
48
            self._energy = self._energy.at(pos)
49 50 51 52 53 54 55

    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
56
        return _toFlatNdarray(self._energy.gradient)
57 58 59

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


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.
    """

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

    def __call__(self, energy):
        import scipy.optimize as opt
86 87 88 89 90 91 92 93
        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
94
            else:
95 96 97 98 99 100
                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
101
        if not r.success:
102
            logger.error("Problem in Scipy minimization: {}".format(r.message))
103 104
            return hlp._energy, IterationController.ERROR
        return hlp._energy, IterationController.CONVERGED
Martin Reinecke's avatar
Martin Reinecke committed
105 106


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

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


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

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


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
149
                return _toNdarray(self._op(_toField(inp, self._op.domain)))
150 151

        op = energy._A
Martin Reinecke's avatar
Martin Reinecke committed
152 153
        b = _toNdarray(energy._b)
        sx = _toNdarray(energy.position)
154 155 156 157 158 159 160 161 162 163
        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
164
        return energy.at(_toField(res, op.domain)), stat