nifty_tools.py 43.8 KB
Newer Older
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
## NIFTY (Numerical Information Field Theory) has been developed at the
## Max-Planck-Institute for Astrophysics.
##
## Copyright (C) 2013 Max-Planck-Society
##
## Author: Marco Selig
## Project homepage: <http://www.mpa-garching.mpg.de/ift/nifty/>
##
## 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/>.

"""
    ..                     __   ____   __
    ..                   /__/ /   _/ /  /_
    ..         __ ___    __  /  /_  /   _/  __   __
    ..       /   _   | /  / /   _/ /  /   /  / /  /
    ..      /  / /  / /  / /  /   /  /_  /  /_/  /
    ..     /__/ /__/ /__/ /__/    \___/  \___   /  tools
    ..                                  /______/

Marco Selig's avatar
Marco Selig committed
31
32
33
34
35
36
    This module extends NIFTY with a nifty set of tools including further
    operators, namely the :py:class:`invertible_operator` and the
    :py:class:`propagator_operator`, and minimization schemes, namely
    :py:class:`steepest_descent` and :py:class:`conjugate_gradient`. Those
    tools are supposed to support the user in solving information field
    theoretical problems (almost) without numerical pain.
37
38
39
40

"""
from __future__ import division
#import numpy as np
41
from nifty_core import *
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56


##-----------------------------------------------------------------------------

class invertible_operator(operator):
    """
        ..       __                                       __     __   __        __
        ..     /__/                                     /  /_  /__/ /  /      /  /
        ..     __   __ ___  __   __   _______   _____  /   _/  __  /  /___   /  /   _______
        ..   /  / /   _   ||  |/  / /   __  / /   __/ /  /   /  / /   _   | /  /  /   __  /
        ..  /  / /  / /  / |     / /  /____/ /  /    /  /_  /  / /  /_/  / /  /_ /  /____/
        .. /__/ /__/ /__/  |____/  \______/ /__/     \___/ /__/  \______/  \___/ \______/  operator class

        NIFTY subclass for invertible, self-adjoint (linear) operators

Marco Selig's avatar
Marco Selig committed
57
58
59
60
        The invertible operator class is an abstract class for self-adjoint or
        symmetric (linear) operators from which other more specific operator
        subclassescan be derived. Such operators inherit an automated inversion
        routine, namely conjugate gradient.
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82

        Parameters
        ----------
        domain : space
            The space wherein valid arguments live.
        uni : bool, *optional*
            Indicates whether the operator is unitary or not.
            (default: False)
        imp : bool, *optional*
            Indicates whether the incorporation of volume weights in
            multiplications is already implemented in the `multiply`
            instance methods or not (default: False).
        para : {single object, tuple/list of objects}, *optional*
            This is a freeform tuple/list of parameters that derivatives of
            the operator class can use (default: None).

        See Also
        --------
        operator

        Notes
        -----
Marco Selig's avatar
Marco Selig committed
83
84
85
        This class is not meant to be instantiated. Operator classes derived
        from this one only need a `_multiply` or `_inverse_multiply` instance
        method to perform the other. However, one of them needs to be defined.
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

        Attributes
        ----------
        domain : space
            The space wherein valid arguments live.
        sym : bool
            Indicates whether the operator is self-adjoint or not.
        uni : bool
            Indicates whether the operator is unitary or not.
        imp : bool
            Indicates whether the incorporation of volume weights in
            multiplications is already implemented in the `multiply`
            instance methods or not.
        target : space
            The space wherein the operator output lives.
        para : {single object, list of objects}
            This is a freeform tuple/list of parameters that derivatives of
            the operator class can use. Not used in the base operators.

    """
    def __init__(self,domain,uni=False,imp=False,para=None):
        """
            Sets the standard operator properties.

            Parameters
            ----------
            domain : space
                The space wherein valid arguments live.
            uni : bool, *optional*
                Indicates whether the operator is unitary or not.
                (default: False)
            imp : bool, *optional*
                Indicates whether the incorporation of volume weights in
                multiplications is already implemented in the `multiply`
                instance methods or not (default: False).
            para : {single object, tuple/list of objects}, *optional*
                This is a freeform tuple/list of parameters that derivatives of
                the operator class can use (default: None).

        """
        if(not isinstance(domain,space)):
            raise TypeError(about._errors.cstring("ERROR: invalid input."))
        self.domain = domain
        self.sym = True
        self.uni = bool(uni)

        if(self.domain.discrete):
            self.imp = True
        else:
            self.imp = bool(imp)

137
        self.target = self.domain
138
139
140
141
142
143

        if(para is not None):
            self.para = para

    ##+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

144
    def _multiply(self,x,force=False,W=None,spam=None,reset=None,note=False,x0=None,tol=1E-4,clevel=1,limii=None,**kwargs):
145
        """
146
147
            Applies the invertible operator to a given field by invoking a
            conjugate gradient.
148
149
150

            Parameters
            ----------
151
152
            x : field
                Valid input field.
153
154
155
156
157
158
            force : bool
                Indicates wheter to return a field instead of ``None`` is
                forced incase the conjugate gradient fails.

            Returns
            -------
159
            OIIx : field
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
                Mapped field with suitable domain.

            See Also
            --------
            conjugate_gradient

            Other Parameters
            ----------------
            W : {operator, function}, *optional*
                Operator `W` that is a preconditioner on `A` and is applicable to a
                field (default: None).
            spam : function, *optional*
                Callback function which is given the current `x` and iteration
                counter each iteration (default: None).
            reset : integer, *optional*
                Number of iterations after which to restart; i.e., forget previous
                conjugated directions (default: sqrt(b.dim())).
            note : bool, *optional*
                Indicates whether notes are printed or not (default: False).
            x0 : field, *optional*
                Starting guess for the minimization.
            tol : scalar, *optional*
                Tolerance specifying convergence; measured by current relative
                residual (default: 1E-4).
            clevel : integer, *optional*
                Number of times the tolerance should be undershot before
                exiting (default: 1).
            limii : integer, *optional*
                Maximum number of iterations performed (default: 10 * b.dim()).

        """
191
192
        x_,convergence = conjugate_gradient(self.inverse_times,x,W=W,spam=spam,reset=reset,note=note)(x0=x0,tol=tol,clevel=clevel,limii=limii)
        ## check convergence
193
        if(not convergence):
194
            if(not force)or(x_ is None):
195
196
                return None
            about.warnings.cprint("WARNING: conjugate gradient failed.")
197
198
199
        ## weight if ...
        if(not self.imp): ## continiuos domain/target
            x_.weight(power=-1,overwrite=True)
200
        return x_
201

202
    def _inverse_multiply(self,x,force=False,W=None,spam=None,reset=None,note=False,x0=None,tol=1E-4,clevel=1,limii=None,**kwargs):
203
        """
204
205
            Applies the inverse of the invertible operator to a given field by
            invoking a conjugate gradient.
206
207
208

            Parameters
            ----------
209
210
            x : field
                Valid input field.
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
            force : bool
                Indicates wheter to return a field instead of ``None`` is
                forced incase the conjugate gradient fails.

            Returns
            -------
            OIx : field
                Mapped field with suitable domain.

            See Also
            --------
            conjugate_gradient

            Other Parameters
            ----------------
            W : {operator, function}, *optional*
                Operator `W` that is a preconditioner on `A` and is applicable to a
                field (default: None).
            spam : function, *optional*
                Callback function which is given the current `x` and iteration
                counter each iteration (default: None).
            reset : integer, *optional*
                Number of iterations after which to restart; i.e., forget previous
                conjugated directions (default: sqrt(b.dim())).
            note : bool, *optional*
                Indicates whether notes are printed or not (default: False).
            x0 : field, *optional*
                Starting guess for the minimization.
            tol : scalar, *optional*
                Tolerance specifying convergence; measured by current relative
                residual (default: 1E-4).
            clevel : integer, *optional*
                Number of times the tolerance should be undershot before
                exiting (default: 1).
            limii : integer, *optional*
                Maximum number of iterations performed (default: 10 * b.dim()).

        """
249
250
        x_,convergence = conjugate_gradient(self.times,x,W=W,spam=spam,reset=reset,note=note)(x0=x0,tol=tol,clevel=clevel,limii=limii)
        ## check convergence
251
        if(not convergence):
252
            if(not force)or(x_ is None):
253
254
                return None
            about.warnings.cprint("WARNING: conjugate gradient failed.")
255
256
257
        ## weight if ...
        if(not self.imp): ## continiuos domain/target
            x_.weight(power=1,overwrite=True)
258
        return x_
259

260
261
262
    ##+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

    def __repr__(self):
263
        return "<nifty_tools.invertible_operator>"
264

265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
##-----------------------------------------------------------------------------

##-----------------------------------------------------------------------------

class propagator_operator(operator):
    """
        ..                                                                            __
        ..                                                                          /  /_
        ..      _______   _____   ______    ______    ____ __   ____ __   ____ __  /   _/  ______    _____
        ..    /   _   / /   __/ /   _   | /   _   | /   _   / /   _   / /   _   / /  /   /   _   | /   __/
        ..   /  /_/  / /  /    /  /_/  / /  /_/  / /  /_/  / /  /_/  / /  /_/  / /  /_  /  /_/  / /  /
        ..  /   ____/ /__/     \______/ /   ____/  \______|  \___   /  \______|  \___/  \______/ /__/     operator class
        .. /__/                        /__/                 /______/

        NIFTY subclass for propagator operators (of a certain family)

        The propagator operators :math:`D` implemented here have an inverse
Marco Selig's avatar
Marco Selig committed
282
283
        formulation like :math:`(S^{-1} + M)`, :math:`(S^{-1} + N^{-1})`, or
        :math:`(S^{-1} + R^\dagger N^{-1} R)` as appearing in Wiener filter
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
        theory.

        Parameters
        ----------
        S : operator
            Covariance of the signal prior.
        M : operator
            Likelihood contribution.
        R : operator
            Response operator translating signal to (noiseless) data.
        N : operator
            Covariance of the noise prior or the likelihood, respectively.

        See Also
        --------
        conjugate_gradient

        Notes
        -----
Marco Selig's avatar
Marco Selig committed
303
        The propagator will puzzle the operators `S` and `M` or `R`, `N` or
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
        only `N` together in the predefined from, a domain is set
        automatically. The application of the inverse is done by invoking a
        conjugate gradient.
        Note that changes to `S`, `M`, `R` or `N` auto-update the propagator.

        Examples
        --------
        >>> f = field(rg_space(4), val=[2, 4, 6, 8])
        >>> S = power_operator(f.target, spec=1)
        >>> N = diagonal_operator(f.domain, diag=1)
        >>> D = propagator_operator(S=S, N=N) # D^{-1} = S^{-1} + N^{-1}
        >>> D(f).val
        array([ 1.,  2.,  3.,  4.])

        Attributes
        ----------
        domain : space
            A space wherein valid arguments live.
        codomain : space
            An alternative space wherein valid arguments live; commonly the
            codomain of the `domain` attribute.
        sym : bool
            Indicates that the operator is self-adjoint.
        uni : bool
            Indicates that the operator is not unitary.
        imp : bool
            Indicates that volume weights are implemented in the `multiply`
            instance methods.
        target : space
            The space wherein the operator output lives.
        _A1 : {operator, function}
            Application of :math:`S^{-1}` to a field.
        _A2 : {operator, function}
            Application of all operations not included in `A1` to a field.
        RN : {2-tuple of operators}, *optional*
            Contains `R` and `N` if given.

    """
    def __init__(self,S=None,M=None,R=None,N=None):
        """
            Sets the standard operator properties and `codomain`, `_A1`, `_A2`,
            and `RN` if required.

            Parameters
            ----------
            S : operator
                Covariance of the signal prior.
            M : operator
                Likelihood contribution.
            R : operator
                Response operator translating signal to (noiseless) data.
            N : operator
                Covariance of the noise prior or the likelihood, respectively.

        """
        ## check signal prior covariance
        if(S is None):
            raise Exception(about._errors.cstring("ERROR: insufficient input."))
        elif(not isinstance(S,operator)):
            raise ValueError(about._errors.cstring("ERROR: invalid input."))
        space1 = S.domain

        ## check likelihood (pseudo) covariance
        if(M is None):
            if(N is None):
                raise Exception(about._errors.cstring("ERROR: insufficient input."))
            elif(not isinstance(N,operator)):
                raise ValueError(about._errors.cstring("ERROR: invalid input."))
            if(R is None):
                space2 = N.domain
            elif(not isinstance(R,operator)):
                raise ValueError(about._errors.cstring("ERROR: invalid input."))
            else:
                space2 = R.domain
        elif(not isinstance(M,operator)):
            raise ValueError(about._errors.cstring("ERROR: invalid input."))
        else:
            space2 = M.domain

        ## set spaces
        self.domain = space2
        if(self.domain.check_codomain(space1)):
            self.codomain = space1
        else:
            raise ValueError(about._errors.cstring("ERROR: invalid input."))
        self.target = self.domain

        ## define A1 == S_inverse
        if(isinstance(S,diagonal_operator)):
            self._A1 = S._inverse_multiply ## S.imp == True
        else:
            self._A1 = S.inverse_times

        ## define A2 == M == R_adjoint N_inverse R == N_inverse
        if(M is None):
            if(R is not None):
                self.RN = (R,N)
                if(isinstance(N,diagonal_operator)):
                    self._A2 = self._standard_M_times_1
                else:
                    self._A2 = self._standard_M_times_2
            elif(isinstance(N,diagonal_operator)):
                self._A2 = N._inverse_multiply ## N.imp == True
            else:
                self._A2 = N.inverse_times
        elif(isinstance(M,diagonal_operator)):
            self._A2 = M._multiply ## M.imp == True
        else:
            self._A2 = M.times

        self.sym = True
        self.uni = False
        self.imp = True

    ##+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

420
    def _standard_M_times_1(self,x,**kwargs): ## applies > R_adjoint N_inverse R assuming N is diagonal
421
422
        return self.RN[0].adjoint_times(self.RN[1]._inverse_multiply(self.RN[0].times(x))) ## N.imp = True

423
    def _standard_M_times_2(self,x,**kwargs): ## applies > R_adjoint N_inverse R
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
        return self.RN[0].adjoint_times(self.RN[1].inverse_times(self.RN[0].times(x)))

    ##+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

    def _inverse_multiply_1(self,x,**kwargs): ## > applies A1 + A2 in self.codomain
        return self._A1(x)+self._A2(x.transform(self.domain)).transform(self.codomain)

    def _inverse_multiply_2(self,x,**kwargs): ## > applies A1 + A2 in self.domain
        return self._A1(x.transform(self.codomain)).transform(self.domain)+self._A2(x)

    ##+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

    def _briefing(self,x): ## > prepares x for `multiply`
        ## inspect x
        if(not isinstance(x,field)):
            return field(self.domain,val=x,target=None),False
        ## check x.domain
        elif(x.domain==self.domain):
            return x,False
        elif(x.domain==self.codomain):
            return x,True
        ## transform
        else:
            return x.transform(target=self.codomain,overwrite=False),True

    def _debriefing(self,x,x_,in_codomain): ## > evaluates x and x_ after `multiply`
        if(x_ is None):
            return None
        ## inspect x
        elif(isinstance(x,field)):
            ## repair ...
            if(in_codomain)and(x.domain!=self.codomain):
                    x_ = x_.transform(target=x.domain,overwrite=False) ## ... domain
            if(x_.target!=x.target):
                x_.set_target(newtarget=x.target) ## ... codomain
        return x_

    ##+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

    def times(self,x,force=False,W=None,spam=None,reset=None,note=False,x0=None,tol=1E-4,clevel=1,limii=None,**kwargs):
        """
465
466
            Applies the propagator to a given object by invoking a
            conjugate gradient.
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517

            Parameters
            ----------
            x : {scalar, list, array, field}
                Scalars are interpreted as constant arrays, and an array will
                be interpreted as a field on the domain of the operator.
            force : bool
                Indicates wheter to return a field instead of ``None`` is
                forced incase the conjugate gradient fails.

            Returns
            -------
            Dx : field
                Mapped field with suitable domain.

            See Also
            --------
            conjugate_gradient

            Other Parameters
            ----------------
            W : {operator, function}, *optional*
                Operator `W` that is a preconditioner on `A` and is applicable to a
                field (default: None).
            spam : function, *optional*
                Callback function which is given the current `x` and iteration
                counter each iteration (default: None).
            reset : integer, *optional*
                Number of iterations after which to restart; i.e., forget previous
                conjugated directions (default: sqrt(b.dim())).
            note : bool, *optional*
                Indicates whether notes are printed or not (default: False).
            x0 : field, *optional*
                Starting guess for the minimization.
            tol : scalar, *optional*
                Tolerance specifying convergence; measured by current relative
                residual (default: 1E-4).
            clevel : integer, *optional*
                Number of times the tolerance should be undershot before
                exiting (default: 1).
            limii : integer, *optional*
                Maximum number of iterations performed (default: 10 * b.dim()).

        """
        ## prepare
        x_,in_codomain = self._briefing(x)
        ## apply operator
        if(in_codomain):
            A = self._inverse_multiply_1
        else:
            A = self._inverse_multiply_2
518
        x_,convergence = conjugate_gradient(A,x_,W=W,spam=spam,reset=reset,note=note)(x0=x0,tol=tol,clevel=clevel,limii=limii)
519
        ## evaluate
520
        if(not convergence):
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
            if(not force):
                return None
            about.warnings.cprint("WARNING: conjugate gradient failed.")
        return self._debriefing(x,x_,in_codomain)

    def inverse_times(self,x,**kwargs):
        """
            Applies the inverse propagator to a given object.

            Parameters
            ----------
            x : {scalar, list, array, field}
                Scalars are interpreted as constant arrays, and an array will
                be interpreted as a field on the domain of the operator.

            Returns
            -------
            DIx : field
                Mapped field with suitable domain.

        """
        ## prepare
        x_,in_codomain = self._briefing(x)
        ## apply operator
        if(in_codomain):
            x_ = self._inverse_multiply_1(x_)
        else:
            x_ = self._inverse_multiply_2(x_)
        ## evaluate
        return self._debriefing(x,x_,in_codomain)

552
553
554
    ##+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

    def __repr__(self):
555
        return "<nifty_tools.propagator_operator>"
556

557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
##-----------------------------------------------------------------------------

##=============================================================================

class conjugate_gradient(object):
    """
        ..      _______       ____ __
        ..    /  _____/     /   _   /
        ..   /  /____  __  /  /_/  / __
        ..   \______//__/  \____  //__/  class
        ..                /______/

        NIFTY tool class for conjugate gradient

        This tool minimizes :math:`A x = b` with respect to `x` given `A` and
        `b` using a conjugate gradient; i.e., a step-by-step minimization
        relying on conjugated gradient directions. Further, `A` is assumed to
        be a positive definite and self-adjoint operator. The use of a
        preconditioner `W` that is roughly the inverse of `A` is optional.
        For details on the methodology refer to [#]_, for details on usage and
        output, see the notes below.

        Parameters
        ----------
        A : {operator, function}
            Operator `A` applicable to a field.
        b : field
            Resulting field of the operation `A(x)`.
        W : {operator, function}, *optional*
            Operator `W` that is a preconditioner on `A` and is applicable to a
            field (default: None).
        spam : function, *optional*
            Callback function which is given the current `x` and iteration
            counter each iteration (default: None).
        reset : integer, *optional*
            Number of iterations after which to restart; i.e., forget previous
            conjugated directions (default: sqrt(b.dim())).
        note : bool, *optional*
            Indicates whether notes are printed or not (default: False).

597
598
599
600
        See Also
        --------
        scipy.sparse.linalg.cg

601
602
603
604
605
606
607
        Notes
        -----
        After initialization by `__init__`, the minimizer is started by calling
        it using `__call__`, which takes additional parameters. Notifications,
        if enabled, will state the iteration number, current step widths
        `alpha` and `beta`, the current relative residual `delta` that is
        compared to the tolerance, and the convergence level if changed.
608
609
610
611
612
        The minimizer will exit in three states: DEAD if alpha becomes
        infinite, QUIT if the maximum number of iterations is reached, or DONE
        if convergence is achieved. Returned will be the latest `x` and the
        latest convergence level, which can evaluate ``True`` for the exit
        states QUIT and DONE.
613
614
615
616

        References
        ----------
        .. [#] J. R. Shewchuk, 1994, `"An Introduction to the Conjugate
Marco Selig's avatar
Marco Selig committed
617
618
            Gradient Method Without the Agonizing Pain"
            <http://www.cs.cmu.edu/~quake-papers/painless-conjugate-gradient.pdf>`_
619
620
621
622
623

        Examples
        --------
        >>> b = field(point_space(2), val=[1, 9])
        >>> A = diagonal_operator(b.domain, diag=[4, 3])
624
        >>> x,convergence = conjugate_gradient(A, b, note=True)(tol=1E-4, clevel=3)
625
626
627
628
629
630
        iteration : 00000001   alpha = 3.3E-01   beta = 1.3E-03   delta = 3.6E-02
        iteration : 00000002   alpha = 2.5E-01   beta = 7.6E-04   delta = 1.0E-03
        iteration : 00000003   alpha = 3.3E-01   beta = 2.5E-04   delta = 1.6E-05   convergence level : 1
        iteration : 00000004   alpha = 2.5E-01   beta = 1.8E-06   delta = 2.1E-08   convergence level : 2
        iteration : 00000005   alpha = 2.5E-01   beta = 2.2E-03   delta = 1.0E-09   convergence level : 3
        ... done.
631
        >>> bool(convergence)
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
        True
        >>> x.val # yields 1/4 and 9/3
        array([ 0.25,  3.  ])

        Attributes
        ----------
        A : {operator, function}
            Operator `A` applicable to a field.
        x : field
            Current field.
        b : field
            Resulting field of the operation `A(x)`.
        W : {operator, function}
            Operator `W` that is a preconditioner on `A` and is applicable to a
            field; can be ``None``.
        spam : function
            Callback function which is given the current `x` and iteration
            counter each iteration; can be ``None``.
        reset : integer
            Number of iterations after which to restart; i.e., forget previous
            conjugated directions (default: sqrt(b.dim())).
        note : notification
            Notification instance.

    """
    def __init__(self,A,b,W=None,spam=None,reset=None,note=False):
        """
            Initializes the conjugate_gradient and sets the attributes (except
            for `x`).

            Parameters
            ----------
            A : {operator, function}
                Operator `A` applicable to a field.
            b : field
                Resulting field of the operation `A(x)`.
            W : {operator, function}, *optional*
                Operator `W` that is a preconditioner on `A` and is applicable to a
                field (default: None).
            spam : function, *optional*
                Callback function which is given the current `x` and iteration
                counter each iteration (default: None).
            reset : integer, *optional*
                Number of iterations after which to restart; i.e., forget previous
                conjugated directions (default: sqrt(b.dim())).
            note : bool, *optional*
                Indicates whether notes are printed or not (default: False).

        """
        if(hasattr(A,"__call__")):
            self.A = A ## applies A
        else:
            raise AttributeError(about._errors.cstring("ERROR: invalid input."))
        self.b = b

        if(W is None)or(hasattr(W,"__call__")):
            self.W = W ## applies W ~ A_inverse
        else:
            raise AttributeError(about._errors.cstring("ERROR: invalid input."))

        self.spam = spam ## serves as callback given x and iteration number
        if(reset is None): ## 2 < reset ~ sqrt(dim)
            self.reset = max(2,int(np.sqrt(b.domain.dim(split=False))))
        else:
            self.reset = max(2,int(reset))
        self.note = notification(default=bool(note))

    ##+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

    def __call__(self,x0=None,**kwargs): ## > runs cg with/without preconditioner
        """
            Runs the conjugate gradient minimization.

            Parameters
            ----------
            x0 : field, *optional*
                Starting guess for the minimization.
            tol : scalar, *optional*
                Tolerance specifying convergence; measured by current relative
                residual (default: 1E-4).
            clevel : integer, *optional*
                Number of times the tolerance should be undershot before
                exiting (default: 1).
            limii : integer, *optional*
                Maximum number of iterations performed (default: 10 * b.dim()).

            Returns
            -------
            x : field
                Latest `x` of the minimization.
722
723
724
            convergence : integer
                Latest convergence level indicating whether the minimization
                has converged or not.
725
726
727
728
729
730
731
732
733
734
735
736

        """
        self.x = field(self.b.domain,val=x0,target=self.b.target)

        if(self.W is None):
            return self._calc_without(**kwargs)
        else:
            return self._calc_with(**kwargs)

    ##+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

    def _calc_without(self,tol=1E-4,clevel=1,limii=None): ## > runs cg without preconditioner
737
        clevel = int(clevel)
738
739
        if(limii is None):
            limii = 10*self.b.domain.dim(split=False)
740
741
        else:
            limii = int(limii)
742

743
744
        r = self.b-self.A(self.x)
        d = field(self.b.domain,val=np.copy(r.val),target=self.b.target)
745
        gamma = r.dot(d)
746
747
        if(gamma==0):
            return self.x,clevel+1
748
749
750
751
752
753
754
        delta_ = np.absolute(gamma)**(-0.5)

        convergence = 0
        ii = 1
        while(True):
            q = self.A(d)
            alpha = gamma/d.dot(q) ## positive definite
755
756
757
            if(not np.isfinite(alpha)):
                self.note.cprint("\niteration : %08u   alpha = NAN\n... dead."%ii)
                return self.x,0
758
            self.x += alpha*d
759
760
761
762
            if(np.signbit(np.real(alpha))):
                about.warnings.cprint("WARNING: positive definiteness of A violated.")
                r = self.b-self.A(self.x)
            elif(ii%self.reset==0):
763
764
765
766
767
768
769
770
771
                r = self.b-self.A(self.x)
            else:
                r -= alpha*q
            gamma_ = gamma
            gamma = r.dot(r)
            beta = max(0,gamma/gamma_) ## positive definite
            d = r+beta*d

            delta = delta_*np.absolute(gamma)**0.5
Marco Selig's avatar
Marco Selig committed
772
            self.note.cflush("\niteration : %08u   alpha = %3.1E   beta = %3.1E   delta = %3.1E"%(ii,np.real(alpha),np.real(beta),np.real(delta)))
773
774
775
            if(ii==limii):
                self.note.cprint("\n... quit.")
                break
776
            elif(gamma==0):
777
                convergence = clevel+1
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
                self.note.cprint("   convergence level : INF\n... done.")
                break
            elif(np.absolute(delta)<tol):
                convergence += 1
                self.note.cflush("   convergence level : %u"%convergence)
                if(convergence==clevel):
                    self.note.cprint("\n... done.")
                    break
            else:
                convergence = max(0,convergence-1)

            if(self.spam is not None):
                self.spam(self.x,ii)

            ii += 1

        if(self.spam is not None):
            self.spam(self.x,ii)

797
        return self.x,convergence
798
799
800
801
802

    ##+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

    def _calc_with(self,tol=1E-4,clevel=1,limii=None): ## > runs cg with preconditioner

803
        clevel = int(clevel)
804
805
        if(limii is None):
            limii = 10*self.b.domain.dim(split=False)
806
807
        else:
            limii = int(limii)
808
809
810
811

        r = self.b-self.A(self.x)
        d = self.W(r)
        gamma = r.dot(d)
812
813
        if(gamma==0):
            return self.x,clevel+1
814
815
816
817
818
819
820
        delta_ = np.absolute(gamma)**(-0.5)

        convergence = 0
        ii = 1
        while(True):
            q = self.A(d)
            alpha = gamma/d.dot(q) ## positive definite
821
822
823
            if(not np.isfinite(alpha)):
                self.note.cprint("\niteration : %08u   alpha = NAN\n... dead."%ii)
                return self.x,0
824
            self.x += alpha*d ## update
825
826
827
828
            if(np.signbit(np.real(alpha))):
                about.warnings.cprint("WARNING: positive definiteness of A violated.")
                r = self.b-self.A(self.x)
            elif(ii%self.reset==0):
829
830
831
832
833
834
835
836
837
838
                r = self.b-self.A(self.x)
            else:
                r -= alpha*q
            s = self.W(r)
            gamma_ = gamma
            gamma = r.dot(s)
            beta = max(0,gamma/gamma_) ## positive definite
            d = s+beta*d ## conjugated gradient

            delta = delta_*np.absolute(gamma)**0.5
Marco Selig's avatar
Marco Selig committed
839
            self.note.cflush("\niteration : %08u   alpha = %3.1E   beta = %3.1E   delta = %3.1E"%(ii,np.real(alpha),np.real(beta),np.real(delta)))
840
841
842
            if(ii==limii):
                self.note.cprint("\n... quit.")
                break
843
844
            elif(gamma==0):
                convergence = clevel+1
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
                self.note.cprint("   convergence level : INF\n... done.")
                break
            elif(np.absolute(delta)<tol):
                convergence += 1
                self.note.cflush("   convergence level : %u"%convergence)
                if(convergence==clevel):
                    self.note.cprint("\n... done.")
                    break
            else:
                convergence = max(0,convergence-1)

            if(self.spam is not None):
                self.spam(self.x,ii)

            ii += 1

        if(self.spam is not None):
            self.spam(self.x,ii)

864
        return self.x,convergence
865

866
867
868
    ##+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

    def __repr__(self):
869
        return "<nifty_tools.conjugate_gradient>"
870

871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
##=============================================================================





##=============================================================================

class steepest_descent(object):
    """
        ..                          __
        ..                        /  /
        ..      _______      ____/  /
        ..    /  _____/    /   _   /
        ..   /_____  / __ /  /_/  / __
        ..  /_______//__/ \______|/__/  class

        NIFTY tool class for steepest descent minimization

        This tool minimizes a scalar energy-function by steepest descent using
        the functions gradient. Steps and step widths are choosen according to
        the Wolfe conditions [#]_. For details on usage and output, see the
        notes below.

        Parameters
        ----------
        eggs : function
            Given the current `x` it returns the tuple of energy and gradient.
        spam : function, *optional*
            Callback function which is given the current `x` and iteration
            counter each iteration (default: None).
        a : {4-tuple}, *optional*
            Numbers obeying 0 < a1 ~ a2 < 1 ~ a3 < a4 that modify the step
            widths (default: (0.2,0.5,1,2)).
        c : {2-tuple}, *optional*
            Numbers obeying 0 < c1 < c2 < 1 that specify the Wolfe-conditions
            (default: (1E-4,0.9)).
        note : bool, *optional*
            Indicates whether notes are printed or not (default: False).

911
912
913
914
915
        See Also
        --------
        scipy.optimize.fmin_cg, scipy.optimize.fmin_ncg,
        scipy.optimize.fmin_l_bfgs_b

916
917
918
919
920
921
922
923
924
        Notes
        -----
        After initialization by `__init__`, the minimizer is started by calling
        it using `__call__`, which takes additional parameters. Notifications,
        if enabled, will state the iteration number, current step width `alpha`,
        current maximal change `delta` that is compared to the tolerance, and
        the convergence level if changed. The minimizer will exit in three
        states: DEAD if no step width above 1E-13 is accepted, QUIT if the
        maximum number of iterations is reached, or DONE if convergence is
925
926
        achieved. Returned will be the latest `x` and the latest convergence
        level, which can evaluate ``True`` for all exit states.
927
928
929
930
931
932
933
934
935
936
937
938
939
940

        References
        ----------
        .. [#] J. Nocedal and S. J. Wright, Springer 2006, "Numerical
            Optimization", ISBN: 978-0-387-30303-1 (print) / 978-0-387-40065-5
            `(online) <http://link.springer.com/book/10.1007/978-0-387-40065-5/page/1>`_

        Examples
        --------
        >>> def egg(x):
        ...     E = 0.5*x.dot(x) # energy E(x) -- a two-dimensional parabola
        ...     g = x # gradient
        ...     return E,g
        >>> x = field(point_space(2), val=[1, 3])
941
        >>> x,convergence = steepest_descent(egg, note=True)(x0=x, tol=1E-4, clevel=3)
942
943
944
945
946
947
948
949
        iteration : 00000001   alpha = 1.0E+00   delta = 6.5E-01
        iteration : 00000002   alpha = 2.0E+00   delta = 1.4E-01
        iteration : 00000003   alpha = 1.6E-01   delta = 2.1E-03
        iteration : 00000004   alpha = 2.6E-03   delta = 3.0E-04
        iteration : 00000005   alpha = 2.0E-04   delta = 5.3E-05   convergence level : 1
        iteration : 00000006   alpha = 8.2E-05   delta = 4.4E-06   convergence level : 2
        iteration : 00000007   alpha = 6.6E-06   delta = 3.1E-06   convergence level : 3
        ... done.
950
        >>> bool(convergence)
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
        True
        >>> x.val # approximately zero
        array([ -6.87299426e-07  -2.06189828e-06])

        Attributes
        ----------
        x : field
            Current field.
        eggs : function
            Given the current `x` it returns the tuple of energy and gradient.
        spam : function
            Callback function which is given the current `x` and iteration
            counter each iteration; can be ``None``.
        a : {4-tuple}
            Numbers obeying 0 < a1 ~ a2 < 1 ~ a3 < a4 that modify the step
            widths (default: (0.2,0.5,1,2)).
        c : {2-tuple}
            Numbers obeying 0 < c1 < c2 < 1 that specify the Wolfe-conditions
            (default: (1E-4,0.9)).
        note : notification
            Notification instance.

    """
    def __init__(self,eggs,spam=None,a=(0.2,0.5,1,2),c=(1E-4,0.9),note=False):
        """
            Initializes the steepest_descent and sets the attributes (except
            for `x`).

            Parameters
            ----------
            eggs : function
                Given the current `x` it returns the tuple of energy and gradient.
            spam : function, *optional*
                Callback function which is given the current `x` and iteration
                counter each iteration (default: None).
            a : {4-tuple}, *optional*
                Numbers obeying 0 < a1 ~ a2 < 1 ~ a3 < a4 that modify the step
                widths (default: (0.2,0.5,1,2)).
            c : {2-tuple}, *optional*
                Numbers obeying 0 < c1 < c2 < 1 that specify the Wolfe-conditions
                (default: (1E-4,0.9)).
            note : bool, *optional*
                Indicates whether notes are printed or not (default: False).

        """
        self.eggs = eggs ## returns energy and gradient

        self.spam = spam ## serves as callback given x and iteration number
        self.a = a ## 0 < a1 ~ a2 < 1 ~ a3 < a4
        self.c = c ## 0 < c1 < c2 < 1
        self.note = notification(default=bool(note))

    ##+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

    def __call__(self,x0,alpha=1,tol=1E-4,clevel=8,limii=100000):
        """
            Runs the steepest descent minimization.

            Parameters
            ----------
            x0 : field
                Starting guess for the minimization.
            alpha : scalar, *optional*
                Starting step width to be multiplied with normalized gradient
                (default: 1).
            tol : scalar, *optional*
                Tolerance specifying convergence; measured by maximal change in
                `x` (default: 1E-4).
            clevel : integer, *optional*
                Number of times the tolerance should be undershot before
                exiting (default: 8).
            limii : integer, *optional*
                Maximum number of iterations performed (default: 100,000).

            Returns
            -------
            x : field
                Latest `x` of the minimization.
1029
1030
1031
            convergence : integer
                Latest convergence level indicating whether the minimization
                has converged or not.
1032
1033
1034
1035
1036

        """
        if(not isinstance(x0,field)):
            raise TypeError(about._errors.cstring("ERROR: invalid input."))
        self.x = x0
1037
1038
1039
1040

        clevel = int(clevel)
        limii = int(limii)

1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
        E,g = self.eggs(self.x) ## energy and gradient
        norm = g.norm() ## gradient norm

        convergence = 0
        ii = 1
        while(True):
            x_,E,g,alpha,a = self._get_alpha(E,g,norm,alpha) ## "news",alpha,a

            if(alpha is None):
                self.note.cprint("\niteration : %08u   alpha < 1.0E-13\n... dead."%ii)
                break
            else:
                delta = np.absolute(g.val).max()*(alpha/norm)
                self.note.cflush("\niteration : %08u   alpha = %3.1E   delta = %3.1E"%(ii,alpha,delta))
                ## update
                self.x = x_
                alpha *= a

            norm = g.norm() ## gradient norm
            if(ii==limii):
                self.note.cprint("\n... quit.")
                break
            elif(delta<tol):
                convergence += 1
                self.note.cflush("   convergence level : %u"%convergence)
                if(convergence==clevel):
1067
                    convergence += int(ii==clevel)
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
                    self.note.cprint("\n... done.")
                    break
            else:
                convergence = max(0,convergence-1)

            if(self.spam is not None):
                self.spam(self.x,ii)

            ii += 1

        if(self.spam is not None):
            self.spam(self.x,ii)

1081
        return self.x,convergence
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129

    ##+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

    def _get_alpha(self,E,g,norm,alpha): ## > determines the new alpha
        while(True):
            ## Wolfe conditions
            wolfe,x_,E_,g_,a = self._check_wolfe(E,g,norm,alpha)
#            wolfe,x_,E_,g_,a = self._check_strong_wolfe(E,g,norm,alpha)
            if(wolfe):
                return x_,E_,g_,alpha,a
            else:
                alpha *= a
                if(alpha<1E-13):
                    return None,None,None,None,None

    ##+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

    def _check_wolfe(self,E,g,norm,alpha): ## > checks the Wolfe conditions
        x_ = self._get_x(g,norm,alpha)
        pg = norm
        E_,g_ = self.eggs(x_)
        if(E_>E+self.c[0]*alpha*pg):
            if(E_<E):
                return True,x_,E_,g_,self.a[1]
            return False,None,None,None,self.a[0]
        pg_ = g.dot(g_)/norm
        if(pg_<self.c[1]*pg):
            return True,x_,E_,g_,self.a[3]
        return True,x_,E_,g_,self.a[2]

#    def _check_strong_wolfe(self,E,g,norm,alpha): ## > checks the strong Wolfe conditions
#        x_ = self._get_x(g,norm,alpha)
#        pg = norm
#        E_,g_ = self.eggs(x_)
#        if(E_>E+self.c[0]*alpha*pg):
#            if(E_<E):
#                return True,x_,E_,g_,self.a[1]
#            return False,None,None,None,self.a[0]
#        apg_ = np.absolute(g.dot(g_))/norm
#        if(apg_>self.c[1]*np.absolute(pg)):
#            return True,x_,E_,g_,self.a[3]
#        return True,x_,E_,g_,self.a[2]

    ##+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

    def _get_x(self,g,norm,alpha): ## > updates x
        return self.x-g*(alpha/norm)

1130
1131
1132
    ##+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++

    def __repr__(self):
1133
        return "<nifty_tools.steepest_descent>"
1134

1135
1136
##=============================================================================