demo_excaliwir.py 10.6 KB
Newer Older
Marco Selig's avatar
Marco Selig committed
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
## 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/>.

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

    NIFTY demo for (critical) Wiener filtering of Gaussian random signals.

"""
from __future__ import division
from nifty import *
from scipy.sparse.linalg import LinearOperator as lo
from scipy.sparse.linalg import cg


note = notification()


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

## spaces
r1 = rg_space(512,1,zerocenter=False)
r2 = rg_space(64,2)
h = hp_space(16)
g = gl_space(48)
z = s_space = k = k_space = p = d_space = None

## power spectrum (and more)
53
power = kindex = rho = powerindex = powerundex = None
Marco Selig's avatar
Marco Selig committed
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

## operators
S = Sk = R = N = Nj = D = None

## fields
s = n = d = j = m = None

## propagator class
class propagator_operator(operator):
    """
        This is the information propagator from the Wiener filter formula.
        It is defined by its inverse. It is given the prior signal covariance S,
        the noise covariance N and the response R in para.
    """
    def _inverse_multiply(self,x):
        ## The inverse can be calculated directly
        S,N,R = self.para
        return S.inverse_times(x)+R.adjoint_times(N.inverse_times(R.times(x)))

    ## the inverse multiplication and multiplication with S modified to return 1d arrays
    _matvec = (lambda self,x: self.inverse_times(x).val.flatten())
    _precon = (lambda self,x: self.para[0].times(x).val.flatten())

    def _multiply(self,x):
        """
            the operator is defined by its inverse, so multiplication has to be
            done by inverting the inverse numerically using the conjugate gradient
            method from scipy
        """
        A = lo(shape=tuple(self.dim()),matvec=self._matvec,dtype=self.domain.datatype) ## linear operator
        b = x.val.flatten()
        x_,info = cg(A,b,x0=None,tol=1.0E-5,maxiter=10*len(x),xtype=None,M=None,callback=None) ## conjugate gradient
        if(info==0):
            return x_
        else:
            note.cprint("NOTE: conjugate gradient failed.")
            return None

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



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

def setup(space,s2n=3,nvar=None):
    """
        sets up the spaces, operators and fields

        Parameters
        ----------
        space : space
            the signal lives in `space`
        s2n : positive number, *optional*
            `s2n` is the signal to noise ratio (default: 3)
        nvar = positive number, *optional*
            the noise variance, `nvar` will be calculated according to
            `s2n` if not specified (default: None)
    """
112
    global z,s_space,k,k_space,p,d_space,power,kindex,rho,powerindex,powerundex,S,Sk,R,N,Nj,D,s,n,d,j,m
Marco Selig's avatar
Marco Selig committed
113
114
115
116
117
118
119

    ## signal space
    z = s_space = space

    ## conjugate space
    k = k_space = s_space.get_codomain()
    ## the power indices are calculated once and saved
120
    kindex,rho,powerindex,powerundex = k_space.get_power_indices()
Marco Selig's avatar
Marco Selig committed
121
122

    ## power spectrum
123
    power = (lambda kk: 42/(kk+1)**3)
Marco Selig's avatar
Marco Selig committed
124
125

    ## prior signal covariance operator (power operator)
126
    S = power_operator(k_space,spec=power)
Marco Selig's avatar
Marco Selig committed
127
    ## projection operator to the spectral bands
128
    Sk = S.get_projection_operator()
Marco Selig's avatar
Marco Selig committed
129
    ## the Gaussian random field generated from its power operator S
130
    s = S.get_random_field(domain=s_space,target=k_space)
Marco Selig's avatar
Marco Selig committed
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

    ## response
    R = response_operator(s_space,sigma=0,mask=1)
    ## data space
    p = d_space = R.target

    ## calculating the noise covariance
    if(nvar is None):
        svar = np.var(s.val) ## given unit response
        nvar = svar/s2n**2
    ## noise covariance operator
    N = diagonal_operator(d_space,diag=nvar,bare=True)
    ## Gaussian noise generated from its covariance N
    n = N.get_random_field(domain=d_space,target=d_space)

    ## data
    d = R(s)+n

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

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

def run(space=r1,s2n=3,nvar=None,**kwargs):
    """
        runs the demo of the generalised Wiener filter

        Parameters
        ----------
        space : space, *optional*
            `space` can be any space from nifty, that supports the plotting
            routine (default: r1 = rg_space(512,1,zerocenter=False))
        s2n : positive number, *optional*
            `s2n` is the signal to noise (default: 3)
        nvar : positive number, *optional*
            the noise variance, `nvar` will be calculated according to
            `s2n` if not specified (default: None)
    """
168
    global s_space,k_space,d_space,power,S,Sk,R,N,Nj,D,s,n,d,j,m
Marco Selig's avatar
Marco Selig committed
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

    ## setting up signal, noise, data and the operators S, N and R
    setup(space,s2n=s2n,nvar=nvar)

    ## information source
    j = R.adjoint_times(N.inverse_times(d))

    ## information propagator
    D = propagator_operator(s_space,sym=True,imp=True,para=[S,N,R])

    ## reconstructed map
    m = D(j)
    if(m is None):
        return None

    ## fields
    s.plot(title="signal",**kwargs)
#    n.cast_domain(s_space,newtarget=k_space)
#    n.plot(title="noise",**kwargs)
#    n.cast_domain(d_space,newtarget=d_space)
    d.cast_domain(s_space,newtarget=k_space)
    d.plot(title="data",vmin=np.min(s.val),vmax=np.max(s.val),**kwargs)
    d.cast_domain(d_space,newtarget=d_space)
    m.plot(title="reconstructed map",vmin=np.min(s.val),vmax=np.max(s.val),**kwargs)

    ## power spectrum
195
#    s.plot(title="power spectra",power=True,other=(m,power),mono=False)
Marco Selig's avatar
Marco Selig committed
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226

    ## uncertainty
#    uncertainty = D.hat(bare=True,nrun=D.domain.dim()//4,target=k_space)
#    if(np.all(uncertainty.val>0)):
#        sqrt(uncertainty).plot(title="standard deviation",**kwargs)

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

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

def run_critical(space=r2,s2n=3,nvar=None,q=1E-12,alpha=1,perception=[1,0],**kwargs):
    """
        runs the demo of the critical generalised Wiener filter

        Parameters
        ----------
        space : space, *optional*
            `space` can be any space from nifty, that supports the plotting
            routine (default: r2 = rg_space(64,2))
        s2n : positive number, *optional*
            `s2n` is the signal to noise (default: 3)
        nvar : positive number, *optional*
            the noise variance, `nvar` will be calculated according to
            `s2n` if not specified (default: None)
        q : positive number, *optional*
            `q` is the minimal power on all scales (default: 1E-12)
        alpha : a number >= 1, *optional*
            `alpha` = 1 means Jeffreys prior for the power spectrum (default: 1)
        perception : array of shape (2,1), *optional*
            perception[0] is delta, perception[1] is epsilon. They are tuning
            factors for the filter (default: [1,0])
227
228
229
230
231

        See Also
        --------
        infer_power

Marco Selig's avatar
Marco Selig committed
232
    """
233
    global s_space,k_space,d_space,power,rho,S,Sk,R,N,Nj,D,s,n,d,j,m
Marco Selig's avatar
Marco Selig committed
234
235
236
237
238
239
240
241
242
243

    ## setting up signal, noise, data and the operators S, N and R
    setup(space,s2n=s2n,nvar=nvar)
    if(perception[1] is None):
        perception[1] = rho/2*(perception[0]-1)

    ## information source
    j = R.adjoint_times(N.inverse_times(d))

    ## unknown power spectrum, the power operator is given an initial value
244
    S.set_power(42,bare=True) ## The answer is 42. I double checked.
Marco Selig's avatar
Marco Selig committed
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
    ## the power spectrum is drawn from the first guess power operator
    pk = S.get_power(bare=False) ## non-bare(!)

    ## information propagator
    D = propagator_operator(s_space,sym=True,imp=True,para=[S,N,R])

    ## iterative reconstruction of the power spectrum and the map
    iteration = 0
    while(True):
        iteration += 1

        ## the Wiener filter reconstruction using the current power spectrum
        m = D(j)
        if(m is None):
            return None

        ## measuring a new estimated power spectrum from the current reconstruction
        b1 = Sk.pseudo_tr(m) ## == Sk(m).pseudo_dot(m), but faster
        b2 = Sk.pseudo_tr(D,nrun=np.sqrt(Sk.domain.dim())//4) ## probing of the partial traces of D
        pk_new = (2*q+b1+perception[0]*b2)/(rho+2*(alpha-1+perception[1])) ## non-bare(!)
265
        pk_new = smooth_power(pk_new,domain=k_space,mode="2s",exclude=min(8,len(rho))) ## smoothing
Marco Selig's avatar
Marco Selig committed
266
        ## the power operator is given the new spectrum
267
        S.set_power(pk_new,bare=False) ## auto-updates D
Marco Selig's avatar
Marco Selig committed
268
269
270
271

        ## check convergence
        log_change = np.max(np.abs(log(pk_new/pk)))
        if(log_change<=0.01):
272
            note.cprint("NOTE: desired accuracy reached in iteration %u."%(iteration))
Marco Selig's avatar
Marco Selig committed
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
            break
        else:
            note.cprint("NOTE: log-change == %4.3f ( > 1%% ) in iteration %u."%(log_change,iteration))
            pk = np.copy(pk_new)

    ## fields
    s.plot(title="signal",**kwargs)
#    n.cast_domain(s_space,newtarget=k_space)
#    n.plot(title="noise",**kwargs)
#    n.cast_domain(d_space,newtarget=d_space)
    d.cast_domain(s_space,newtarget=k_space)
    d.plot(title="data",vmin=np.min(s.val),vmax=np.max(s.val),**kwargs)
    d.cast_domain(d_space,newtarget=d_space)
    m.plot(title="reconstructed map",vmin=np.min(s.val),vmax=np.max(s.val),**kwargs)

    ## power spectrum
289
    s.plot(title="power spectra",power=True,other=(S.get_power(),power),mono=False)
Marco Selig's avatar
Marco Selig committed
290
291
292
293
294
295
296
297

    ## uncertainty
#    uncertainty = D.hat(bare=True,nrun=D.domain.dim()//4,target=k_space)
#    if(np.all(uncertainty.val>0)):
#        sqrt(uncertainty).plot(title="standard deviation")

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