inversion.py 45.6 KB
Newer Older
Maximilian Schanner's avatar
Maximilian Schanner committed
1
"""
Maximilian Schanner's avatar
Maximilian Schanner committed
2
3
    This module is part of the CORBASS algorithm. It can be seen as the core
    module, where the actual inversion code is defined.
Maximilian Schanner's avatar
Maximilian Schanner committed
4

Maximilian Schanner's avatar
Maximilian Schanner committed
5
6
    Copyright (C) 2019 Helmholtz Centre Potsdam GFZ,
        German Research Centre for Geosciences, Potsdam, Germany
Maximilian Schanner's avatar
Maximilian Schanner committed
7
8

    Cite as:
Maximilian Schanner's avatar
Maximilian Schanner committed
9
    Schanner, Maximilian Arthus and Mauerberger, Stefan (2019)
Maximilian Schanner's avatar
Maximilian Schanner committed
10
    CORBASS: CORrelation Based Archeomagnetic SnapShot model. V. 1.0.
Maximilian Schanner's avatar
Maximilian Schanner committed
11
    GFZ Data Services. http://doi.org/10.5880/GFZ.2.3.2019.008
Maximilian Schanner's avatar
Maximilian Schanner committed
12

Maximilian Schanner's avatar
Maximilian Schanner committed
13
14
    This file is part of CORBASS.

Maximilian Schanner's avatar
Maximilian Schanner committed
15
16
17
18
19
20
21
22
23
24
25
    CORBASS 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.

    CORBASS 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
Maximilian Schanner's avatar
Maximilian Schanner committed
26
    along with CORBASS. If not, see <https://www.gnu.org/licenses/>.
Maximilian Schanner's avatar
Maximilian Schanner committed
27
"""
Maximilian Schanner's avatar
Maximilian Schanner committed
28
29
30
if __name__ == '__main__':
    # make relative imports work
    import sys
Maximilian Schanner's avatar
Maximilian Schanner committed
31
32
    from pathlib import Path
    sys.path.append(str(Path(__file__).absolute().parents[1]))
Maximilian Schanner's avatar
Maximilian Schanner committed
33
34
35
36
37
38

import numpy as np
from pandas import DataFrame

import pyfield

39
from corbass.utils import scaling, read_data, dif2nez, nez2dif
Maximilian Schanner's avatar
Maximilian Schanner committed
40

Maximilian Schanner's avatar
Maximilian Schanner committed
41
RCMB = 3480.
Maximilian Schanner's avatar
Maximilian Schanner committed
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


def grad_d(b):
    """ Calculate the gradient of D """
    res = np.column_stack((-b[1::3], b[0::3], np.zeros_like(b[0::3])))
    res /= (b[0::3]**2 + b[1::3]**2)[:, np.newaxis]
    return np.rad2deg(res.flatten())


def grad_i(b):
    """ Calculate the gradient of I """
    b = np.asarray(b)
    res = -np.column_stack((b[0::3], b[1::3], b[2::3])) \
        * b[2::3][:, np.newaxis]
    res /= (b[0::3]**2 + b[1::3]**2 + b[2::3]**2)[:, np.newaxis]
    res[:, 2] += 1.
    res /= np.sqrt(b[0::3]**2 + b[1::3]**2)[:, np.newaxis]
    return np.rad2deg(res.flatten())


def grad_f(b):
    """ Calculate the gradient of F """
    b = np.asarray(b)
    res = np.column_stack((b[0::3], b[1::3], b[2::3]))
    res /= np.sqrt(b[0::3]**2 + b[1::3]**2 + b[2::3]**2)[:, np.newaxis]
    return res.flatten()


class Inversion:
71
    def __init__(self, fname, drop_duplicates=False):
Maximilian Schanner's avatar
Maximilian Schanner committed
72
73
74
75
76
77
        """
        Parameters
        ----------
        fname : str or pandas.DataFrame
            Either a file path (str) leading to or a pandas.DataFrame
            containing archeomagnetic data
78
79
80
        drop_duplicates : bool (optional, default is False)
            Wether to drop multiple records for the same location and keep only
            the first one.
Maximilian Schanner's avatar
Maximilian Schanner committed
81
82
83
        """
        self._get_obs_locs(fname, drop_duplicates=drop_duplicates)

84
    def full_run(self, r_Ref, lamb, epsilon, rho, to_calc):
Maximilian Schanner's avatar
Maximilian Schanner committed
85
86
87
88
89
90
91
92
        """ Calculate quantities specified in to_calc using a single call of
        _prep_calc. The idea is to save calculation time since likelihood and
        inversion share several calculation steps.

        Parameters
        ----------
        r_Ref : float
            The reference radius.
93
        lamb : float
Maximilian Schanner's avatar
Maximilian Schanner committed
94
            The scaling factor for the non-dipole kernel contribution.
95
        epsilon : float
Maximilian Schanner's avatar
Maximilian Schanner committed
96
            Scaling factor for the observational error.
97
        rho : float
Maximilian Schanner's avatar
Maximilian Schanner committed
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
            A residual error.
        to_calc : dict
            A dictionary whose keys specify which quantities to calculate and
            whose items give the corresponding parameters. Possible keys are:
            'field' : n_desi (int)
                Invert for the field, parameter is the number of design
                points.
            'dif' : n_desi (int)
                Invert for DIF, parameter is the number of design points.
                If both 'field' and 'dif' are in the dict, n_desi from 'field'
                is used.
            'coeffs' : lmax (int)
                Invert for the Gauss coefficients, parameter is the highest
                degree.
            'lik' : diag (bool)
                Calculate the likelihood, parameter is wether to return the
                individual terms for diagnosis and plotting.
115
116
117
            'zf' : n_desi (int)
                Invert for the down component and intensity at the CMB,
                parameter is the number of design points.
Maximilian Schanner's avatar
Maximilian Schanner committed
118
119
120
121
122
123
124
125
126
127
128
            'pst' : prior (callable) or None
                Calculate the posterior for the kernel-parameters, parameter is
                either a callable of the signature pri(float, float, float),
                returning the log prior value for the parameters or None, in
                which case the standard prior will be used.

        Returns
        -------
        res : dict
            A dictionary of results. The keys are the same as in to_calc,
            with suffix '_mean' to access the posterior mean and '_cov' to
129
            access the covariances (except for 'pst').
Maximilian Schanner's avatar
Maximilian Schanner committed
130
        """
131
        self._prep_calc(r_Ref, lamb, epsilon, rho)
Maximilian Schanner's avatar
Maximilian Schanner committed
132
133
134
135

        res = {}

        if 'field' in to_calc.keys():
136
            mean, cov = self._field_inversion(lamb, to_calc['field'])
Maximilian Schanner's avatar
Maximilian Schanner committed
137
138
139
140
141
142
            res['field_mean'] = mean
            res['field_cov'] = cov
            if 'dif' in to_calc.keys():
                self.mu_NEZ_desi = mean
                self.cov_NEZ_desi = cov
        elif 'dif' in to_calc.keys():
143
            mean, cov = self._field_inversion(lamb, to_calc['dif'])
Maximilian Schanner's avatar
Maximilian Schanner committed
144
145
146
            self.mu_NEZ_desi = mean
            self.cov_NEZ_desi = cov
        if 'dif' in to_calc.keys():
Maximilian Schanner's avatar
Maximilian Schanner committed
147
148
149
            d, i, f = self._dif_inversion()
            res['dif_mean'] = (d[0], i[0], f[0])
            res['dif_cov'] = (d[1], i[1], f[1])
Maximilian Schanner's avatar
Maximilian Schanner committed
150
        if 'coeffs' in to_calc.keys():
151
            mean, cov = self._coeff_inversion(r_Ref, lamb,
Maximilian Schanner's avatar
Maximilian Schanner committed
152
153
154
155
156
157
158
159
160
161
162
163
164
                                              to_calc['coeffs'])

            ls = [pyfield.i2lm_l(i) for i in range(mean.size)]
            ms = [pyfield.i2lm_m(i) for i in range(mean.size)]

            res['coeffs_mean'] = mean
            res['coeffs_cov'] = cov
            res['coeffs_ls'] = ls
            res['coeffs_ms'] = ms
        if 'lik' in to_calc.keys():
            res['lik'] = self._log_likeli(to_calc['lik'])
        if 'pst' in to_calc.keys():
            if to_calc['pst'] is not None:
165
                res['pst'] = self._log_pst(lamb, epsilon, rho,
Maximilian Schanner's avatar
Maximilian Schanner committed
166
167
                                           log_pri=to_calc['pst'])
            else:
168
                res['pst'] = self._log_pst(lamb, epsilon, rho)
Maximilian Schanner's avatar
Maximilian Schanner committed
169
        if 'zf' in to_calc.keys():
170
            dummy = self._zf_CMB_inversion(lamb, to_calc['zf'])
Maximilian Schanner's avatar
Maximilian Schanner committed
171
172
173
174
175
176
            res['zf_loc'] = dummy[0]
            res['zf_mean'] = (dummy[1][0], dummy[2][0])
            res['zf_var'] = (dummy[1][1], dummy[2][1])

        return res

177
    def _prep_calc(self, r_Ref, lamb, epsilon, rho):
Maximilian Schanner's avatar
Maximilian Schanner committed
178
179
180
181
182
183
184
185
        """ This function serves as a setup for inversion and likelihood.
        During a run several quantities that are used in both functions are
        calculated and set as class attributes.

        Parameters
        ----------
        r_Ref : float
            The reference radius.
186
        lamb : float
Maximilian Schanner's avatar
Maximilian Schanner committed
187
            The scaling factor for the non-dipole kernel contribution.
188
        epsilon : float
Maximilian Schanner's avatar
Maximilian Schanner committed
189
            Scaling factor for the observational error.
190
        rho : float
Maximilian Schanner's avatar
Maximilian Schanner committed
191
192
            A residual error.
        """
193
        # Kernels from the pyfield library
Maximilian Schanner's avatar
Maximilian Schanner committed
194
        self.sp_Dip = pyfield.StatDipIntSpline(pyfield.REARTH)
Maximilian Schanner's avatar
Maximilian Schanner committed
195
196
197
198
        self.sp_WoDip = pyfield.StatWODipIntSpline(r_Ref)
        # Allocate and retrieve basis functions (-grad Phi)
        self.dspharm_obs = np.empty((3, 3*self.n_obs), order='F')
        self.sp_Dip.gaussField(1, self.x_obs, self.dspharm_obs)
199

Maximilian Schanner's avatar
Maximilian Schanner committed
200
201
202
203
204
205
206
        # Prior covariance at locations of observation)
        # Non-dipole part (without dipole)
        cov_NEZ_obs_WoDip = np.zeros((3*self.n_obs, 3*self.n_obs),
                                     order='F')
        self.sp_WoDip.field(self.x_obs, cov_NEZ_obs_WoDip)

        # Construct entire data covariance matrix
207
        self.cov_NEZ_obs = lamb**2*cov_NEZ_obs_WoDip
Maximilian Schanner's avatar
Maximilian Schanner committed
208
209
210
211
212
213
214
215

        # Slice complete records
        dspharm_full = self.dspharm_obs[:, self.idx_NEZ_full]

        self.cov_NEZ_full = self.cov_NEZ_obs[self.idx_NEZ_full[:, None],
                                             self.idx_NEZ_full[None, :]]

        # Add errors for complete records
216
217
        errs = epsilon**2 * self.errs_full \
            + rho**2 * self.epsilon_term_full
Maximilian Schanner's avatar
Maximilian Schanner committed
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

        errs_inv = np.linalg.inv(errs)

        lapl_appr = np.linalg.inv(np.linalg.multi_dot((self.trans_grad.T,
                                                       errs_inv,
                                                       self.trans_grad)))

        self.cov_NEZ_full += lapl_appr

        # Invert covariance using Cholesky
        self.L = np.linalg.cholesky(self.cov_NEZ_full)
        self.L_inv = np.linalg.inv(self.L)

        self.prc_NEZ_full = np.dot(self.L_inv.T, self.L_inv)
        # Non informative stuff...
        self.h = np.dot(self.L_inv, dspharm_full.T)

        self.c = np.dot(self.L_inv, self.o_NEZ_full)

        self.var_g_bar = np.linalg.inv(np.dot(self.h.T, self.h))
        self.g_bar = np.linalg.multi_dot((self.var_g_bar, self.h.T, self.c))

        # Calculate some posteriors that are needed multiple times
        # 1st step
        self.cor_obs_full = np.copy(self.cov_NEZ_obs[:, self.idx_NEZ_full])

        self.mu_NEZ_full = np.dot(self.g_bar,
245
                                  dspharm_full)
Maximilian Schanner's avatar
Maximilian Schanner committed
246
247
248
249
250
251
252
253
254
255
256

        self.mu_NEZ_obs = np.dot(self.g_bar, self.dspharm_obs) \
            + np.linalg.multi_dot([self.cor_obs_full,
                                   self.prc_NEZ_full,
                                   self.o_NEZ_full - self.mu_NEZ_full])

        self.cov_NEZ_obs -= np.linalg.multi_dot([self.cor_obs_full,
                                                 self.prc_NEZ_full,
                                                 self.cor_obs_full.T])

        self.R_obs = self.dspharm_obs \
257
            - np.linalg.multi_dot((dspharm_full,
Maximilian Schanner's avatar
Maximilian Schanner committed
258
259
260
261
262
263
264
265
266
267
268
269
270
271
                                   self.prc_NEZ_full,
                                   self.cor_obs_full.T))

        self.cov_NEZ_obs += np.linalg.multi_dot((self.R_obs.T,
                                                 self.var_g_bar,
                                                 self.R_obs))

        # 2nd step
        # Calculate linearised mean for D, I and F records
        mu_N_obs = self.mu_NEZ_obs[0::3]  # North
        mu_E_obs = self.mu_NEZ_obs[1::3]  # East
        mu_Z_obs = self.mu_NEZ_obs[2::3]  # Down

        # 0th order term
272
273
274
275
276
277
278
279
280
        mu_D_icmp = nez2dif(n=mu_N_obs[self.idx_D],
                            e=mu_E_obs[self.idx_D],
                            d=mu_Z_obs[self.idx_D])[0]
        mu_I_icmp = nez2dif(n=mu_N_obs[self.idx_I],
                            e=mu_E_obs[self.idx_I],
                            d=mu_Z_obs[self.idx_I])[1]
        mu_F_icmp = nez2dif(n=mu_N_obs[self.idx_F],
                            e=mu_E_obs[self.idx_F],
                            d=mu_Z_obs[self.idx_F])[2]
Maximilian Schanner's avatar
Maximilian Schanner committed
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302

        self.mu_DIF_icmp = np.concatenate((mu_D_icmp, mu_I_icmp, mu_F_icmp))

        # 1st order term
        grad_D = grad_d(self.mu_NEZ_obs[self.idx_NEZ_D]).reshape((-1, 3))
        grad_I = grad_i(self.mu_NEZ_obs[self.idx_NEZ_I]).reshape((-1, 3))
        grad_F = grad_f(self.mu_NEZ_obs[self.idx_NEZ_F]).reshape((-1, 3))
        self.grad_DIF = np.concatenate((grad_D, grad_I, grad_F))

        # Covariance w.r.t. north, east and down
        cov_NEZ_icmp = self.cov_NEZ_obs[self.idx_NEZ_icmp[:, None],
                                        self.idx_NEZ_icmp[None, :]]

        cov_NEZ_icmp = cov_NEZ_icmp.reshape(self.n_icmp, 3,
                                            self.n_icmp, 3)

        # so far only BB-covariance, product with grad gives oo-covariance
        self.cov_DIF_icmp = np.einsum('ij,ijlm,lm->il',
                                      self.grad_DIF,
                                      cov_NEZ_icmp,
                                      self.grad_DIF)

303
304
        epsilon_term = np.diag(np.dot(self.grad_DIF, self.grad_DIF.T))
        epsilon_term = np.diag(epsilon_term)
Maximilian Schanner's avatar
Maximilian Schanner committed
305

306
307
        self.cov_DIF_icmp += epsilon**2 * self.errs_icmp \
            + rho**2 * epsilon_term
Maximilian Schanner's avatar
Maximilian Schanner committed
308
309
310
311

        # Precision matrix
        self.prc_DIF_icmp = np.linalg.inv(self.cov_DIF_icmp)

312
    def field_inversion(self, r_Ref, lamb, epsilon, rho, n_desi):
Maximilian Schanner's avatar
Maximilian Schanner committed
313
314
315
316
317
318
        """ Invert for the magnetic field at equidistributed design points.

        Parameters
        ----------
        r_Ref : float
            The reference radius.
319
        lamb : float
Maximilian Schanner's avatar
Maximilian Schanner committed
320
            The scaling factor for the non-dipole kernel contribution.
321
        epsilon : float
Maximilian Schanner's avatar
Maximilian Schanner committed
322
            Scaling factor for the observational error.
323
        rho : float
Maximilian Schanner's avatar
Maximilian Schanner committed
324
325
326
327
328
329
330
331
332
333
334
            A residual error.
        n_desi : int
            The number of design points.

        Returns
        -------
        mu_NEZ_desi : array
            The posterior mean of the field at design points.
        cov_NEZ_desi : array
            The posterior covariance of the field at design points.
        """
335
336
        self._prep_calc(r_Ref, lamb, epsilon, rho)
        return self._field_inversion(lamb, n_desi=n_desi)
Maximilian Schanner's avatar
Maximilian Schanner committed
337

338
    def _field_inversion(self, lamb, n_desi=10):
Maximilian Schanner's avatar
Maximilian Schanner committed
339
340
341
342
343
344
        """ Private version of field_inversion. Should be called only after
        _prep_calc has been called with the parameters of interest, since it
        makes use of the stuff calculated there!

        Parameters
        ----------
345
        lamb : float
Maximilian Schanner's avatar
Maximilian Schanner committed
346
347
348
349
350
351
352
353
354
355
356
357
            The scaling factor for the non-dipole kernel contribution.
        n_desi : int (optional, default is 10)
            The number of design points.

        Returns
        -------
        mu_NEZ_desi : array
            The posterior mean of the field at design points.
        cov_NEZ_desi : array
            The posterior covariance of the field at design points.
        """
        # Setup design points
358
359
360
361
362
363
364
365
366
367
368
369
        if isinstance(n_desi, np.ndarray):
            self.x_desi = n_desi
            self.n_desi = self.x_desi.shape[1]
        elif isinstance(n_desi, DataFrame):
            self.n_desi = len(n_desi)
            self.x_desi = np.zeros(4, self.n_desi)
            self.x_desi[0] = n_desi[['co-lat']]
            self.x_desi[1] = n_desi[['lon']]
            self.x_desi[2] = n_desi[['rad']]
            self.x_desi[3] = n_desi[['t']]
        else:
            self.x_desi, self.n_desi = self.desi_points(n_desi)
Maximilian Schanner's avatar
Maximilian Schanner committed
370
371
372
373
374
375
376
377
378
379
380
381

        # Allocate and retrieve basis functions at design points
        dspharm_desi = np.empty((3, 3*self.n_desi), order='F')
        self.sp_Dip.gaussField(1, self.x_desi, dspharm_desi)

        # Prior covariance at design points
        cov_NEZ_desi_WoDip = np.zeros((3*self.n_desi,
                                       3*self.n_desi),
                                      order='F')
        self.sp_WoDip.field(self.x_desi, cov_NEZ_desi_WoDip)

        # Construct design points covariance matrix
382
        cov_NEZ_desi = lamb**2*cov_NEZ_desi_WoDip
Maximilian Schanner's avatar
Maximilian Schanner committed
383
384
385
386
387
388
389
390
391
392

        # Prior correlations amongst design points and observations
        cor_NEZ_desi_NEZ_obs_WoDip = np.zeros((3*self.n_desi,
                                               3*self.n_obs),
                                              order='F')
        self.sp_WoDip.field(self.x_desi,
                            self.x_obs,
                            cor_NEZ_desi_NEZ_obs_WoDip)

        # Construct design-observation covariance matrix
393
        cor_NEZ_desi_NEZ_obs = lamb**2*cor_NEZ_desi_NEZ_obs_WoDip
Maximilian Schanner's avatar
Maximilian Schanner committed
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
420
421
422
423
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

        # ------------------------------------------------------------------- #
        # 1st step: Full records
        # ------------------------------------------------------------------- #

        # Correlations between design points and full records
        cor_desi_full = cor_NEZ_desi_NEZ_obs[:, self.idx_NEZ_full]

        # Posterior mean
        mu_NEZ_desi = np.dot(self.g_bar, dspharm_desi) \
            + np.linalg.multi_dot([cor_desi_full,
                                   self.prc_NEZ_full,
                                   self.o_NEZ_full - self.mu_NEZ_full])

        # Posterior covariance
        cov_NEZ_desi -= np.linalg.multi_dot([cor_desi_full,
                                             self.prc_NEZ_full,
                                             cor_desi_full.T])

        R_desi = dspharm_desi \
            - np.linalg.multi_dot((self.dspharm_obs[:, self.idx_NEZ_full],
                                   self.prc_NEZ_full,
                                   cor_desi_full.T))

        cov_NEZ_desi += np.linalg.multi_dot((R_desi.T, self.var_g_bar, R_desi))

        # Posterior correlation between design points and observations
        cor_NEZ_desi_NEZ_obs -= np.linalg.multi_dot([cor_desi_full,
                                                     self.prc_NEZ_full,
                                                     self.cor_obs_full.T])

        cor_NEZ_desi_NEZ_obs += np.linalg.multi_dot((R_desi.T,
                                                     self.var_g_bar,
                                                     self.R_obs))
        # ------------------------------------------------------------------- #
        # 2nd step: Incomplete records
        # ------------------------------------------------------------------- #
        # Correlations between design points and incomplete records
        cor_NEZ_desi_NEZ_icmp = cor_NEZ_desi_NEZ_obs[:, self.idx_NEZ_icmp]\
            .reshape(3*self.n_desi, self.n_icmp, 3)

        cor_NEZ_desi_DIF_icmp = np.einsum('ijk,jk->ij',
                                          cor_NEZ_desi_NEZ_icmp,
                                          self.grad_DIF)
        # Posterior mean
        mu_NEZ_desi += np.dot(cor_NEZ_desi_DIF_icmp,
                              np.dot(self.prc_DIF_icmp,
                                     self.o_DIF_icmp - self.mu_DIF_icmp))
        # Posterior covariance
        cov_NEZ_desi -= np.linalg.multi_dot((cor_NEZ_desi_DIF_icmp,
                                             self.prc_DIF_icmp,
                                             cor_NEZ_desi_DIF_icmp.T))

        return mu_NEZ_desi, cov_NEZ_desi

449
    def zf_CMB_inversion(self, r_Ref, lamb, epsilon, rho, n_desi):
Maximilian Schanner's avatar
Maximilian Schanner committed
450
451
452
453
454
455
        """ Invert for the magnetic field at equidistributed design points.

        Parameters
        ----------
        r_Ref : float
            The reference radius.
456
        lamb : float
Maximilian Schanner's avatar
Maximilian Schanner committed
457
            The scaling factor for the non-dipole kernel contribution.
458
        epsilon : float
Maximilian Schanner's avatar
Maximilian Schanner committed
459
            Scaling factor for the observational error.
460
        rho : float
Maximilian Schanner's avatar
Maximilian Schanner committed
461
462
463
464
465
466
            A residual error.
        n_desi : int
            The number of design points.

        Returns
        -------
Maximilian Schanner's avatar
Maximilian Schanner committed
467
468
469
470
471
472
473
474
        x_desi : array
            The location array at which the fields have been evaluated
        (Z, Var_Z) : tuple of arrays
            A tuple containing the mean of the down component at the CMB and
            the corresponding pointwise variance
        (F, Var_F) : tuple of arrays
            A tuple containing the mean of the field intensity at the CMB and
            the corresponding pointwise variance
Maximilian Schanner's avatar
Maximilian Schanner committed
475
        """
476
477
        self._prep_calc(r_Ref, lamb, epsilon, rho)
        return self._zf_CMB_inversion(lamb, n_desi=n_desi)
Maximilian Schanner's avatar
Maximilian Schanner committed
478

479
    def _zf_CMB_inversion(self, lamb, n_desi=10):
Maximilian Schanner's avatar
Typo    
Maximilian Schanner committed
480
        """ Private version of zf_CMB_inversion. Should be called only after
Maximilian Schanner's avatar
Maximilian Schanner committed
481
482
483
484
485
        _prep_calc has been called with the parameters of interest, since it
        makes use of the stuff calculated there!

        Parameters
        ----------
486
        lamb : float
Maximilian Schanner's avatar
Maximilian Schanner committed
487
488
489
490
491
492
            The scaling factor for the non-dipole kernel contribution.
        n_desi : int (optional, default is 10)
            The number of design points.

        Returns
        -------
Maximilian Schanner's avatar
Maximilian Schanner committed
493
494
495
496
497
498
499
500
        x_desi : array
            The location array at which the fields have been evaluated
        (Z, Var_Z) : tuple of arrays
            A tuple containing the mean of the down component at the CMB and
            the corresponding pointwise variance
        (F, Var_F) : tuple of arrays
            A tuple containing the mean of the field intensity at the CMB and
            the corresponding pointwise variance
Maximilian Schanner's avatar
Maximilian Schanner committed
501
502
        """
        # Setup design points
Maximilian Schanner's avatar
Maximilian Schanner committed
503
        x_desi, n_desi = self.desi_points(n_desi, rad=RCMB)
Maximilian Schanner's avatar
Maximilian Schanner committed
504
505
506
507
508
509
510
511
512
513
514
515

        # Allocate and retrieve basis functions at design points
        dspharm_desi = np.empty((3, 3*n_desi), order='F')
        self.sp_Dip.gaussField(1, x_desi, dspharm_desi)

        # Prior covariance at design points
        cov_NEZ_desi_WoDip = np.zeros((3*n_desi,
                                       3*n_desi),
                                      order='F')
        self.sp_WoDip.field(x_desi, cov_NEZ_desi_WoDip)

        # Construct design points covariance matrix
516
        cov_NEZ_desi = lamb**2*cov_NEZ_desi_WoDip
Maximilian Schanner's avatar
Maximilian Schanner committed
517
518
519
520
521
522
523
524
525
526

        # Prior correlations amongst design points and observations
        cor_NEZ_desi_NEZ_obs_WoDip = np.zeros((3*n_desi,
                                               3*self.n_obs),
                                              order='F')
        self.sp_WoDip.field(x_desi,
                            self.x_obs,
                            cor_NEZ_desi_NEZ_obs_WoDip)

        # Construct design-observation covariance matrix
527
        cor_NEZ_desi_NEZ_obs = lamb**2*cor_NEZ_desi_NEZ_obs_WoDip
Maximilian Schanner's avatar
Maximilian Schanner committed
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
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

        # ------------------------------------------------------------------- #
        # 1st step: Full records
        # ------------------------------------------------------------------- #

        # Correlations between design points and full records
        cor_desi_full = cor_NEZ_desi_NEZ_obs[:, self.idx_NEZ_full]

        # Posterior mean
        mu_NEZ_desi = np.dot(self.g_bar, dspharm_desi) \
            + np.linalg.multi_dot([cor_desi_full,
                                   self.prc_NEZ_full,
                                   self.o_NEZ_full - self.mu_NEZ_full])

        # Posterior covariance
        cov_NEZ_desi -= np.linalg.multi_dot([cor_desi_full,
                                             self.prc_NEZ_full,
                                             cor_desi_full.T])

        R_desi = dspharm_desi \
            - np.linalg.multi_dot((self.dspharm_obs[:, self.idx_NEZ_full],
                                   self.prc_NEZ_full,
                                   cor_desi_full.T))

        cov_NEZ_desi += np.linalg.multi_dot((R_desi.T, self.var_g_bar, R_desi))

        # Posterior correlation between design points and observations
        cor_NEZ_desi_NEZ_obs -= np.linalg.multi_dot([cor_desi_full,
                                                     self.prc_NEZ_full,
                                                     self.cor_obs_full.T])

        cor_NEZ_desi_NEZ_obs += np.linalg.multi_dot((R_desi.T,
                                                     self.var_g_bar,
                                                     self.R_obs))
        # ------------------------------------------------------------------- #
        # 2nd step: Incomplete records
        # ------------------------------------------------------------------- #
        # Correlations between design points and incomplete records
        cor_NEZ_desi_NEZ_icmp = cor_NEZ_desi_NEZ_obs[:, self.idx_NEZ_icmp]\
            .reshape(3*n_desi, self.n_icmp, 3)

        cor_NEZ_desi_DIF_icmp = np.einsum('ijk,jk->ij',
                                          cor_NEZ_desi_NEZ_icmp,
                                          self.grad_DIF)
        # Posterior mean
        mu_NEZ_desi += np.dot(cor_NEZ_desi_DIF_icmp,
                              np.dot(self.prc_DIF_icmp,
                                     self.o_DIF_icmp - self.mu_DIF_icmp))
        # Posterior covariance
        cov_NEZ_desi -= np.linalg.multi_dot((cor_NEZ_desi_DIF_icmp,
                                             self.prc_DIF_icmp,
                                             cor_NEZ_desi_DIF_icmp.T))

        Z = mu_NEZ_desi[2::3]
        var_Z = cov_NEZ_desi.diagonal()[2::3]

584
        _, _, F = nez2dif(*mu_NEZ_desi.reshape(-1, 3).T)
Maximilian Schanner's avatar
Maximilian Schanner committed
585
586
587
588
589
590
591
592
593
594
595
596

        grad_F_desi = grad_f(mu_NEZ_desi).reshape(-1, 3)

        cov_F = np.einsum('ij,ijkl,kl->ik',
                          grad_F_desi,
                          cov_NEZ_desi.reshape(n_desi, 3,
                                               n_desi, 3),
                          grad_F_desi)
        var_F = cov_F.diagonal()

        return x_desi, (Z, var_Z), (F, var_F)

597
    def coeff_inversion(self, r_Ref, lamb, epsilon, rho, lmax):
Maximilian Schanner's avatar
Maximilian Schanner committed
598
599
600
601
602
603
        """ Invert for the Gauss coefficients.

        Parameters
        ----------
        r_Ref : float
            The reference radius.
604
        lamb : float
Maximilian Schanner's avatar
Maximilian Schanner committed
605
            The scaling factor for the non-dipole kernel contribution.
606
        epsilon : float
Maximilian Schanner's avatar
Maximilian Schanner committed
607
            Scaling factor for the observational error.
608
        rho : float
Maximilian Schanner's avatar
Maximilian Schanner committed
609
610
611
612
613
614
615
616
617
618
619
            A residual error.
        lmax : int
            The highest degree to invert for.

        Returns
        -------
        mu_g : array
            The posterior mean coefficients.
        cov_g : array
            The posterior coefficients-covariance.
        """
620
621
        self._prep_calc(r_Ref, lamb, epsilon, rho)
        return self._coeff_inversion(r_Ref, lamb, lmax=lmax)
Maximilian Schanner's avatar
Maximilian Schanner committed
622

623
    def _coeff_inversion(self, r_Ref, lamb, lmax=14):
Maximilian Schanner's avatar
Maximilian Schanner committed
624
625
626
627
628
629
630
631
        """ Private version of _coeff_inversion. Should be called only after
        _prep_calc has been called with the parameters of interest, since it
        makes use of the stuff calculated there!

        Parameters
        ----------
        r_Ref : float
            The reference radius.
632
        lamb : float
Maximilian Schanner's avatar
Maximilian Schanner committed
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
            The scaling factor for the non-dipole kernel contribution.
        lmax : int (optional, default is 14)
            The highest degree to invert for.

        Returns
        -------
        mu_g : array
            The posterior mean coefficients.
        cov_g : array
            The posterior coefficients-covariance.
        """
        # Total number of coefficients
        self.n_coeff = pyfield.lmax2N(lmax)

        # Prior mean and covariance
        mu_g = np.zeros(self.n_coeff)
649
        var_g = lamb**2*np.ones(self.n_coeff)
Maximilian Schanner's avatar
Maximilian Schanner committed
650
651
652
653
654
655
656
        cov_g = np.diag(var_g)
        cov_g[0:3, 0:3] = np.zeros((3, 3))

        # Correlations amongst Gauss coefficients and observations
        cor_g_obs_WoDip = np.empty((self.n_coeff, 3*self.n_obs), order='F')
        self.sp_WoDip.gaussField(lmax, self.x_obs, cor_g_obs_WoDip)
        # Sum for inversion
657
        cor_g_obs = lamb**2*cor_g_obs_WoDip
Maximilian Schanner's avatar
Maximilian Schanner committed
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

        # correlations amongst coefficients and full records
        cor_g_full = cor_g_obs[:, self.idx_NEZ_full]

        # ------------------------------------------------------------------- #
        # 1st step: Full records
        # ------------------------------------------------------------------- #
        dspharm_full = self.dspharm_obs[:, self.idx_NEZ_full]
        # Posterior mean
        mu_g += np.linalg.multi_dot((cor_g_full,
                                     self.prc_NEZ_full,
                                     self.o_NEZ_full - self.mu_NEZ_full))
        # Posterior covariance
        cov_g -= np.einsum('ik,kl,jl->ij',
                           cor_g_full,
                           self.prc_NEZ_full,
                           cor_g_full,
                           optimize=True)

        cov_g += np.linalg.multi_dot((cor_g_full,
                                      self.prc_NEZ_full,
                                      dspharm_full.T,
                                      self.var_g_bar,
                                      dspharm_full,
                                      self.prc_NEZ_full,
                                      cor_g_full.T))

        cov_g_g_dip = -np.linalg.multi_dot((cor_g_full[3:, :],
                                            self.prc_NEZ_full,
                                            dspharm_full.T,
                                            self.var_g_bar))

        cov_g[:3, 3:] = cov_g_g_dip.T
        cov_g[3:, :3] = cov_g_g_dip

        # Posterior correlation between coefficients and observations
        cor_g_obs -= np.linalg.multi_dot((cor_g_full,
                                          self.prc_NEZ_full,
                                          self.cor_obs_full.T))

        dummy = np.linalg.multi_dot((dspharm_full,
                                     self.prc_NEZ_full,
                                     self.cor_obs_full.T))
        dummy -= self.dspharm_obs

        cor_g_obs += np.linalg.multi_dot((cor_g_full,
                                          self.prc_NEZ_full,
                                          dspharm_full.T,
                                          self.var_g_bar,
                                          dummy))

        # Dipole part from the non-informative prior
Maximilian Schanner's avatar
Maximilian Schanner committed
710
        sss = scaling(pyfield.REARTH, r_Ref, 1)[0]
Maximilian Schanner's avatar
Maximilian Schanner committed
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
        # Scale to R_ref
        mu_g[0:3] = np.multiply(self.g_bar, sss)
        cov_g[0:3, 0:3] = np.multiply(self.var_g_bar, sss**2)
        dummy = self.dspharm_obs \
            - np.linalg.multi_dot((dspharm_full,
                                   self.prc_NEZ_full,
                                   self.cor_obs_full.T))

        cor_g_obs[0:3, :] = sss*np.dot(self.var_g_bar, dummy)

        # ------------------------------------------------------------------- #
        # 2nd step: Incomplete records
        # ------------------------------------------------------------------- #
        # Linearized correlations amongst coefficients and incomplete records
        cor_g_NEZ_icmp = cor_g_obs[:, self.idx_NEZ_icmp]\
            .reshape(self.n_coeff, self.n_icmp, 3)
        cor_g_DIF_icmp = np.einsum('ijk,jk->ij', cor_g_NEZ_icmp, self.grad_DIF)

        # Posterior
        mu_g += np.linalg.multi_dot((cor_g_DIF_icmp,
                                     self.prc_DIF_icmp,
                                     self.o_DIF_icmp - self.mu_DIF_icmp))

        cov_g -= np.linalg.multi_dot((cor_g_DIF_icmp,
                                      self.prc_DIF_icmp,
                                      cor_g_DIF_icmp.T))

        return mu_g, cov_g

740
    def dif_inversion(self, r_Ref, lamb, epsilon, rho, n_desi):
Maximilian Schanner's avatar
Maximilian Schanner committed
741
742
743
744
745
746
        """ Invert for declination, inclination and intensity.

        Parameters
        ----------
        r_Ref : float
            The reference radius.
747
        lamb : float
Maximilian Schanner's avatar
Maximilian Schanner committed
748
            The scaling factor for the non-dipole kernel contribution.
749
        epsilon : float
Maximilian Schanner's avatar
Maximilian Schanner committed
750
            Scaling factor for the observational error.
751
        rho : float
Maximilian Schanner's avatar
Maximilian Schanner committed
752
753
754
755
756
757
            A residual error.
        n_desi : int
            The number of design points.

        Returns
        -------
Maximilian Schanner's avatar
Maximilian Schanner committed
758
759
760
761
762
763
764
765
766
        (mu_D, cov_D) : tuple of arrays
            A tuple containing the posterior mean and the covariance of the
            declination
        (mu_I, cov_I) : tuple of arrays
            A tuple containing the posterior mean and the covariance of the
            inclination
        (mu_F, cov_F) : tuple of arrays
            A tuple containing the posterior mean and the covariance of the
            intensity
Maximilian Schanner's avatar
Maximilian Schanner committed
767
        """
768
        self._prep_calc(r_Ref, lamb, epsilon, rho)
Maximilian Schanner's avatar
Maximilian Schanner committed
769
        self.mu_NEZ_desi, self.cov_NEZ_desi = \
770
            self._field_inversion(lamb, n_desi=n_desi)
Maximilian Schanner's avatar
Maximilian Schanner committed
771
772
773
774
775
776
777
778
779
780
        return self._dif_inversion()

    def _dif_inversion(self):
        """ Private version of dif_inversion. Should be called only after
        _prep_calc has been called and the posterior field has been set as
        class attributes self.mu_NEZ_desi and self.cov_NEZ_desi at the
        corresponding design points!

        Returns
        -------
Maximilian Schanner's avatar
Maximilian Schanner committed
781
782
783
784
785
786
787
788
789
        (mu_D, cov_D) : tuple of arrays
            A tuple containing the posterior mean and the covariance of the
            declination
        (mu_I, cov_I) : tuple of arrays
            A tuple containing the posterior mean and the covariance of the
            inclination
        (mu_F, cov_F) : tuple of arrays
            A tuple containing the posterior mean and the covariance of the
            intensity
Maximilian Schanner's avatar
Maximilian Schanner committed
790
        """
791
        mu_D, mu_I, mu_F = nez2dif(*self.mu_NEZ_desi.reshape(-1, 3).T)
Maximilian Schanner's avatar
Maximilian Schanner committed
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812

        grad_D_desi = grad_d(self.mu_NEZ_desi).reshape(-1, 3)
        grad_I_desi = grad_i(self.mu_NEZ_desi).reshape(-1, 3)
        grad_F_desi = grad_f(self.mu_NEZ_desi).reshape(-1, 3)

        cov_D = np.einsum('ij,ijkl,kl->ik',
                          grad_D_desi,
                          self.cov_NEZ_desi.reshape(self.n_desi, 3,
                                                    self.n_desi, 3),
                          grad_D_desi)
        cov_I = np.einsum('ij,ijkl,kl->ik',
                          grad_I_desi,
                          self.cov_NEZ_desi.reshape(self.n_desi, 3,
                                                    self.n_desi, 3),
                          grad_I_desi)
        cov_F = np.einsum('ij,ijkl,kl->ik',
                          grad_F_desi,
                          self.cov_NEZ_desi.reshape(self.n_desi, 3,
                                                    self.n_desi, 3),
                          grad_F_desi)

Maximilian Schanner's avatar
Maximilian Schanner committed
813
        return (mu_D, cov_D), (mu_I, cov_I), (mu_F, cov_F)
Maximilian Schanner's avatar
Maximilian Schanner committed
814

815
    def log_likeli(self, r_Ref, lamb, epsilon, rho,
Maximilian Schanner's avatar
Maximilian Schanner committed
816
817
818
819
820
821
822
                   diag=False):
        """ The log-likelihood function for the kernel-parameters.

        Parameters
        ----------
        r_Ref : float
            The reference radius.
823
        lamb : float
Maximilian Schanner's avatar
Maximilian Schanner committed
824
            The scaling factor for the non-dipole kernel contribution.
825
        epsilon : float
Maximilian Schanner's avatar
Maximilian Schanner committed
826
            Scaling factor for the observational error.
827
        rho : float
Maximilian Schanner's avatar
Maximilian Schanner committed
828
829
830
831
832
833
834
835
836
837
            A residual error.
        diag : bool (optional, default is False)
            Wether to return the individual terms together with the likelihood,
            for diagnosis and plotting.

        Returns
        -------
        The log-likelihood for given parameters.
        """
        try:
838
            self._prep_calc(r_Ref, lamb, epsilon, rho)
Maximilian Schanner's avatar
Maximilian Schanner committed
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
        except np.linalg.linalg.LinAlgError:
            # if the matrix is not psd, return log(0) (a large negative value)
            if diag:
                return -1e16, -1e16, -1e16
            else:
                return -1e16

        return self._log_likeli(diag)

    def _log_likeli(self, diag=False):
        """ Private version of the log-likelihood function for the
        kernel-parameters. Should be called only after _prep_calc has been
        called with the parameters of interest, since it makes use of the stuff
        calculated there!

        Parameters
        ----------
        diag : bool (optional, default is False)
            Wether to return the individual terms together with the likelihood,
            for diagnosis and plotting.

        Returns
        -------
        The log-likelihood for given parameters.
        """
        dummy = np.linalg.slogdet(np.dot(self.h.T, self.h))
        logdet = np.sum(2*np.log(np.diag(self.L))) \
Maximilian Schanner's avatar
Maximilian Schanner committed
866
            + dummy[1]
Maximilian Schanner's avatar
Maximilian Schanner committed
867
868
869
870
871
872
873
874

        delta = np.sum(self.c**2) - np.linalg.multi_dot((self.c.T,
                                                         self.h,
                                                         self.var_g_bar,
                                                         self.h.T,
                                                         self.c))

        dummy = np.linalg.slogdet(self.cov_DIF_icmp)
Maximilian Schanner's avatar
Maximilian Schanner committed
875
876
877
878
        if 0 < dummy[0]:
            logdet += dummy[1]
        else:
            raise ValueError("Covariance not positive definite.")
Maximilian Schanner's avatar
Maximilian Schanner committed
879
880
881
882
883
884
885
886
887
888
889
890
891

        delta += np.linalg.multi_dot(((self.o_DIF_icmp - self.mu_DIF_icmp).T,
                                      self.prc_DIF_icmp,
                                      (self.o_DIF_icmp - self.mu_DIF_icmp)))

        if diag:
            return - logdet/2. - delta/2. - (self.n_obs-3)/2.*np.log(2*np.pi),\
                - delta/2.,\
                - logdet/2.
        else:
            return - logdet/2. - delta/2. - (self.n_obs-3)/2.*np.log(2*np.pi) \


892
    def log_pst(self, r_Ref, lamb, epsilon, rho,
Maximilian Schanner's avatar
Maximilian Schanner committed
893
894
895
896
897
898
899
                log_pri=None):
        """ Calculate the log-posterior for the kernel-parameters.

        Parameters
        ----------
        r_Ref : float
            The reference radius.
900
        lamb : float
Maximilian Schanner's avatar
Maximilian Schanner committed
901
            The scaling factor for the non-dipole kernel contribution.
902
        epsilon : float
Maximilian Schanner's avatar
Maximilian Schanner committed
903
            Scaling factor for the observational error.
904
        rho : float
Maximilian Schanner's avatar
Maximilian Schanner committed
905
906
            A residual error.
        pri : callable or None
907
            A prior function over the parameters lamb, epsilon, rho.
Maximilian Schanner's avatar
Maximilian Schanner committed
908
            The signature should be pri(a, s_E, s_R) and pri should return a
909
            float. If None (default), the standard log prior will be used.
Maximilian Schanner's avatar
Maximilian Schanner committed
910
911
912
913
914

        Returns
        -------
        The log-posterior value for the kernel-parameters.
        """
915
916
        self._prep_calc(r_Ref, lamb, epsilon, rho)
        return self._log_pst(lamb, epsilon, rho, log_pri=log_pri)
Maximilian Schanner's avatar
Maximilian Schanner committed
917

918
    def _log_pst(self, lamb, epsilon, rho,
Maximilian Schanner's avatar
Maximilian Schanner committed
919
920
921
922
923
924
925
926
927
928
                 log_pri=None):
        """ Private version of the log-posterior function for the
        kernel-parameters. Should be called only after _prep_calc has been
        called with the parameters of interest, since it makes use of the stuff
        calculated there!

        Parameters
        ----------
        r_Ref : float
            The reference radius.
929
        lamb : float
Maximilian Schanner's avatar
Maximilian Schanner committed
930
            The scaling factor for the non-dipole kernel contribution.
931
        epsilon : float
Maximilian Schanner's avatar
Maximilian Schanner committed
932
            Scaling factor for the observational error.
933
        rho : float
Maximilian Schanner's avatar
Maximilian Schanner committed
934
935
            A residual error.
        pri : callable or None
936
            A prior function over the parameters lamb, epsilon, rho.
Maximilian Schanner's avatar
Maximilian Schanner committed
937
938
939
940
941
942
943
944
945
            The signature should be pri(a, s_E, s_R) and pri should return a
            float. If None, the standard log prior will be used.

        Returns
        -------
        The log-posterior value for the kernel-parameters.
        """
        if log_pri is None:
            log_pri = self.__class__.log_pri
946
        return self._log_likeli() + log_pri(lamb, epsilon, rho)
Maximilian Schanner's avatar
Maximilian Schanner committed
947
948

    @staticmethod
949
    def log_pri(lamb, epsilon, rho):
Maximilian Schanner's avatar
Maximilian Schanner committed
950
951
952
953
954
        """ The standard log prior over the kernel parameters to be used in the
        posterior.

        Parameters
        ----------
955
        lamb : float
Maximilian Schanner's avatar
Maximilian Schanner committed
956
            The scaling factor for the non-dipole kernel contribution.
957
        epsilon : float
Maximilian Schanner's avatar
Maximilian Schanner committed
958
            Scaling factor for the observational error.
959
        rho : float
Maximilian Schanner's avatar
Maximilian Schanner committed
960
961
962
963
964
965
966
            A residual error.

        Returns
        -------
        The log prior value for the given parameters.
        """

967
968
969
        return - np.log(lamb) - np.log(rho)

    def _get_obs_locs(self, fname, drop_duplicates=False):
Maximilian Schanner's avatar
Maximilian Schanner committed
970
971
972
973
974
975
976
977
978
        """ Convenience function to get observations and locations from a file,
        setup indices for complete and incomplete records and calculate the
        approximate error matrix for the complete records.

        Parameters
        ----------
        fname : str or pandas.DataFrame
            Either a file path (str) leading to or a pandas.DataFrame
            containing archeomagnetic data
979
980
981
        drop_duplicates : bool (optional, default is False)
            Wether to drop multiple records for the same location and keep only
            the first one
Maximilian Schanner's avatar
Maximilian Schanner committed
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
        """
        if isinstance(fname, DataFrame):
            # read D,I,F-data
            data = fname.copy()
        else:
            data = read_data(fname)

        # Drop duplicates
        if drop_duplicates:
            data.drop_duplicates(subset=['lat', 'lon'], inplace=True,
                                 keep='first')
        # Reset index
        data.reset_index(inplace=True)

        # XXX Map declination to [-180:180] degrees
        #     I doubt that approach being smart
        data['D'].where(data['D'] <= 180, data['D']-360, inplace=True)
        self.n_obs = len(data)
        self.x_obs = np.asfortranarray(data[['co-lat', 'lon', 'rad', 't']]).T

        # DataFrame indices for full records
        self.idx_full = data.dropna(axis=0, subset=['D', 'I', 'F']).index

        # Number of full records
        self.n_full = self.idx_full.size

1008
        # Indices in the order as expected by the pyfield library
Maximilian Schanner's avatar
Maximilian Schanner committed
1009
1010
1011
1012
        self.idx_NEZ_full = np.vstack([3*self.idx_full,
                                       3*self.idx_full+1,
                                       3*self.idx_full+2]).T.ravel()
        # Transform DIF data to NEZ
1013
1014
1015
        o_NEZ_full = dif2nez(d=data['D'].loc[self.idx_full],
                             i=data['I'].loc[self.idx_full],
                             f=data['F'].loc[self.idx_full])
Maximilian Schanner's avatar
Maximilian Schanner committed
1016
1017
1018
1019
1020

        # Full records
        self.o_NEZ_full = np.asfortranarray(o_NEZ_full).T.ravel()

        # Errors for full records by Laplace approximation
1021
1022
        epsilon = np.zeros((3*self.n_full, 3*self.n_full))
        rho = np.zeros((3*self.n_full, 3*self.n_full))
Maximilian Schanner's avatar
Maximilian Schanner committed
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
        trans_grad = np.zeros((3*self.n_full, 3*self.n_full))

        for it in range(self.n_full):
            jac = np.array([grad_d(self.o_NEZ_full[3*it:3*it+3]),
                            grad_i(self.o_NEZ_full[3*it:3*it+3]),
                            grad_f(self.o_NEZ_full[3*it:3*it+3])])

            errs = np.diag([data['dD'].loc[self.idx_full[it]]**2,
                            data['dI'].loc[self.idx_full[it]]**2,
                            data['dF'].loc[self.idx_full[it]]**2])

1034
            epsilons = np.dot(jac, jac.T)
Maximilian Schanner's avatar
Maximilian Schanner committed
1035
1036

            trans_grad[3*it:3*it+3, 3*it:3*it+3] = jac
1037
1038
            epsilon[3*it:3*it+3, 3*it:3*it+3] = errs
            rho[3*it:3*it+3, 3*it:3*it+3] = epsilons
Maximilian Schanner's avatar
Maximilian Schanner committed
1039
1040

        self.trans_grad = trans_grad
1041
1042
        self.errs_full = epsilon
        self.epsilon_term_full = rho
Maximilian Schanner's avatar
Maximilian Schanner committed
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064

        # DataFrame indices for D, I and F records
        self.idx_D = data.query('D==D and not ( D==D and I==I and F==F )') \
            .index
        self.idx_I = data.query('I==I and not ( D==D and I==I and F==F )') \
            .index
        self.idx_F = data.query('F==F and not ( D==D and I==I and F==F )') \
            .index
        # Number of incomplete records
        self.n_icmp = self.idx_D.size + self.idx_I.size + self.idx_F.size

        # Incomplete records
        self.o_DIF_icmp = np.concatenate((data['D'].loc[self.idx_D],
                                          data['I'].loc[self.idx_I],
                                          data['F'].loc[self.idx_F]))

        # Errors for incomplete records
        std_DIF_icmp = np.concatenate((data['dD'].loc[self.idx_D],
                                       data['dI'].loc[self.idx_I],
                                       data['dF'].loc[self.idx_F]))
        self.errs_icmp = np.diag(std_DIF_icmp**2)

1065
        # Index-foo to have index arrays w.r.t the pyfield order
Maximilian Schanner's avatar
Maximilian Schanner committed
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
        self.idx_NEZ_D = np.vstack(((3*self.idx_D,
                                     3*self.idx_D+1,
                                     3*self.idx_D+2))).T.ravel()

        self.idx_NEZ_I = np.vstack(((3*self.idx_I,
                                     3*self.idx_I+1,
                                     3*self.idx_I+2))).T.ravel()

        self.idx_NEZ_F = np.vstack(((3*self.idx_F,
                                     3*self.idx_F+1,
                                     3*self.idx_F+2))).T.ravel()

        # More index foo
        self.idx_NEZ_icmp = np.concatenate((self.idx_NEZ_D,
                                            self.idx_NEZ_I,
                                            self.idx_NEZ_F))

Maximilian Schanner's avatar
Maximilian Schanner committed
1083
    def desi_points(self, n, rad=pyfield.REARTH, ret_n=True):
Maximilian Schanner's avatar
Maximilian Schanner committed
1084
        """ Return equidistant design points on the sphere, as expected by
1085
        the pyfield library.
Maximilian Schanner's avatar
Maximilian Schanner committed
1086
1087
1088
1089
1090

        Parameters
        ----------
        n : int
            Approximate number of design points.
1091
1092
        rad : float (optional, default is pyfield.REARTH)
            The radius at which the points are generated.
Maximilian Schanner's avatar
Maximilian Schanner committed
1093
1094
1095
1096
1097
1098
        ret_n : bool (optional, default is True)
            Wether to return the actual number of design points.

        Returns
        -------
        x_desi : array-like
1099
1100
            An array containing n design points, organized as expected by the
            pyfield library.
Maximilian Schanner's avatar
Maximilian Schanner committed
1101
1102
        """
        theta, phi = zip(*pyfield.equi_sph(n))
Stefan Mauerberger's avatar
Stefan Mauerberger committed
1103
1104
1105
1106
1107
        # Phi's range is [0,2*pi]. To achieve smooth global plots with cartopy
        # the interval depends on the central longitude (cl) chosen. AFAIK it
        # shall be [-180+ct,180+ct] to have points that hit the projection's
        # edges. If we agree on cl=0, subtracting pi shall do the job.
        phi = np.subtract(phi, np.pi)
Maximilian Schanner's avatar
Maximilian Schanner committed
1108
        colat, lon = np.rad2deg((theta, phi))
1109
1110
        colat[0] += 1.e-6       # BUGFIX pyfield returns NaN at the north pole

Maximilian Schanner's avatar
Maximilian Schanner committed
1111
        n_desi = len(colat)  # Actual number of design points
1112
        # Allocate array and store locations as expected by the pyfield library
Maximilian Schanner's avatar
Maximilian Schanner committed
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
        x_desi = np.asfortranarray((colat,
                                    lon,
                                    np.full(n_desi, rad),
                                    np.zeros(n_desi)))
        if ret_n:
            return x_desi, n_desi
        else:
            return x_desi


if __name__ == '__main__':
    from pathlib import Path

    from matplotlib import pyplot as plt
    import cartopy.crs as ccrs

Maximilian Schanner's avatar
Maximilian Schanner committed
1129
    from corbass.utils import load
Maximilian Schanner's avatar
Maximilian Schanner committed
1130
1131
1132
1133
1134

    pars = load(sys.argv)
    # ungroup the data
    data = pars.data.filter(lambda x: True)

Maximilian Schanner's avatar
Maximilian Schanner committed
1135
1136
1137
1138
1139
1140
1141
1142
1143
    # cut the bin from the data
    it_bin = 1550

    cur_dat = data[(it_bin-50 < data['t']) &
                   (data['t'] <= it_bin+50)].copy()
    # set up inversion class
    test = Inversion(cur_dat)
    # run the inversion for the fields

Maximilian Schanner's avatar
Maximilian Schanner committed
1144
    glob_proj = ccrs.Mollweide()
Maximilian Schanner's avatar
Maximilian Schanner committed
1145
1146
1147
1148
1149
1150
    fig, ax = plt.subplots(2, 2,
                           figsize=(13, 5),
                           subplot_kw={'projection':
                                       glob_proj},
                           squeeze=False)

Maximilian Schanner's avatar
Maximilian Schanner committed
1151
    r_ref = 3000
1152
1153
1154
    lamb = 16000
    epsilon = 1.34
    rho = 5000
Maximilian Schanner's avatar
Maximilian Schanner committed
1155

1156
    x_zf, zs, fs = test.zf_CMB_inversion(r_ref, lamb, epsilon, rho,
Maximilian Schanner's avatar
Maximilian Schanner committed
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
                                         1000)
    lat_zf = 90-x_zf[0]
    lon_zf = x_zf[1]
    lat, lon, dummy = glob_proj.transform_points(ccrs.Geodetic(),
                                                 lon_zf,
                                                 lat_zf).T
    z, var_z = zs
    f, var_f = fs

    ax[0, 0].tripcolor(lat, lon, z)
    ax[0, 0].coastlines(zorder=4)
    ax[0, 0].set_title('Z @ CMB')

    ax[1, 0].tripcolor(lat, lon, var_z,
Maximilian Schanner's avatar
Maximilian Schanner committed
1171
                       cmap=plt.cm.binary)
Maximilian Schanner's avatar
Maximilian Schanner committed
1172
1173
1174
1175
1176
1177
1178
    ax[1, 0].coastlines(zorder=4)

    ax[0, 1].tripcolor(lat, lon, f)
    ax[0, 1].coastlines(zorder=4)
    ax[0, 1].set_title('F @ CMB')

    ax[1, 1].tripcolor(lat, lon, var_f,
Maximilian Schanner's avatar
Maximilian Schanner committed
1179
                       cmap=plt.cm.binary)
Maximilian Schanner's avatar
Maximilian Schanner committed
1180
1181
1182
    ax[1, 1].coastlines(zorder=4)

    plt.show()