pynkowski.data.so3data

  1import numpy as np
  2import healpy as hp
  3from .base_da import DataField
  4from .utils_da import get_theta, healpix_derivatives, healpix_second_derivatives
  5from ..stats.minkowski import _MF_prefactor
  6try:
  7    from tqdm.auto import tqdm
  8except:
  9    tqdm = lambda x: x
 10    print('tqdm not loaded')
 11
 12
 13
 14# This function is meant to be used as a callable class
 15# The class is initialized with two arrays, Q and U.
 16# The function __call__ returns the value of f at a given
 17# position psi.
 18class QUarray(np.ndarray):
 19    '''Array to store Q and U maps in the SO(3) formalism.
 20
 21    Parameters
 22    ----------
 23    Q : np.array
 24        Values of the Q field. It can be an array of arrays to represent several maps.
 25        
 26    U : np.array
 27        Values of the U field. Must have the same shape as `Q`.
 28        
 29    Notes
 30    -----
 31    The class is a subclass of `np.ndarray`, so it can be used as a numpy array.
 32    CAREFUL: the multiplication of two SO3arrays is not the same as the multiplication of the $f$ that they represent. 
 33    In order to multiply two SO3arrays, you have to use the __call__ method, i.e. multiply the two SO3arrays after evaluating them at a given positions psi.
 34    however, you can multiply a SO3array with a scalar, or add or substract SO3arrays, and the result will be the correct SO3array
 35    '''
 36    def __new__(cls, Q, U):
 37        obj = np.asarray([Q,U]).view(cls)
 38        return obj
 39
 40    def __array_finalize__(self, obj):
 41        if obj is None: return
 42
 43    def __call__(self, psi):
 44        '''Return the value of f at a given position psi.
 45        
 46        Parameters
 47        ----------
 48        psi : float or np.array
 49            The position at which to evaluate f, in radians. It must be broadcastable with the shape of `Q` and `U`.
 50            
 51        Returns
 52        -------
 53        f : np.array
 54            The value of f at the given position psi.
 55        '''
 56        return np.asarray(self[0]*np.cos(2.*psi) - self[1]*np.sin(2.*psi))
 57    
 58    def derpsi(self):
 59        '''Return the derivative of f with respect to psi, which can be exactly computed in the SO(3) formalism.
 60
 61        Returns
 62        -------
 63        df : np.array
 64            The derivative of f with respect to psi.
 65        '''
 66        return QUarray(-2.*self[1], 2.*self[0])
 67    
 68    def modulus(self):
 69        '''Return the modulus of f.
 70
 71        Returns
 72        -------
 73        fmod : np.array
 74            The modulus of f.
 75        '''
 76        return np.sqrt(self[0]**2 + self[1]**2)
 77    
 78    def pol_angle(self):
 79        '''Return the polarization angle of f.
 80
 81        Returns
 82        -------
 83        pol_angle : np.array
 84            The polarization angle of f.
 85        '''
 86        angle = np.arctan2(-self[1], self[0])/2
 87        angle = angle + np.pi/2*(1-np.sign(self(angle)))/2  # I don't remember why we did this, but it seems to work
 88        angle[angle<0] = angle[angle<0] + np.pi
 89        return angle
 90
 91
 92class SO3DataField(DataField):
 93    """Class for spherical spin fields in the SO(3) formalism, with Q and U in HEALPix format.
 94
 95    Parameters
 96    ----------
 97    Q : np.array
 98        Values of the Q field.
 99    
100    U : np.array
101        Values of the U field. Must have the same shape as `Q`.
102        
103    name : str, optional
104        Name of the field.
105        Defaul : `'SO(3) DataField'`
106    
107    first_der : np.array or None, optional
108        First **covariant** derivatives of the field in an orthonormal basis of the space. Same structure as `field`, and shape `(3, field.shape)`.
109    
110    second_der : np.array or None, optional
111        Second **covariant** derivatives of the field in an orthonormal basis of the space. Same structure as `field`, and shape `(6, field.shape)`.
112        The order of the derivatives is diagonal first, e.g. in `dim=3`: `11`, `22`, `33`, `12`, `13`, `23`.
113        
114    mask : np.array or None, optional
115        Mask where the field is considered. It is a bool array of the same shape that `Q` and `U`.
116        Default: all data is included.
117        
118    Attributes
119    ----------
120    field : np.array
121        Array of shape `(2, Q.shape)` containing the Q and U values. It can be called as `field(psi)` to obtain the value of $f$ for the given `psi`.
122        
123    dim : int
124        Dimension of the space where the field is defined. In this case, the space is SO(3) and this is 3.
125        
126    name : str
127        Name of the field.
128    
129    first_der : list or None
130        First **covariant** derivatives of the field in an orthonormal basis of the space. It a list of size `3`, each entry has the same structure as `field`.
131    
132    second_der : list or None
133        Second **covariant** derivatives of the field in an orthonormal basis of the space. It a list of size `6`, each entry has the same structure as `field`.
134        The order of the derivatives is diagonal first, i.e., `11`, `22`, `33`, `12`, `13`, `23`.
135        
136    mask : np.array
137        Sky mask where the field is considered. It is a bool array of the same shape that `Q` and `U`.
138
139    """   
140    def __init__(self, Q, U, name="SO(3) DataField", first_der=None, second_der=None, mask=None):
141        assert Q.shape == U.shape, "`Q` and `U` must have the same shape"
142        field = QUarray(Q, U)
143        super().__init__(field, 3, name=name, mask=mask, first_der=first_der, second_der=second_der)
144        if mask is None:
145            self.mask = np.ones(Q.shape, dtype=bool)
146        else:
147            assert mask.shape == Q.shape, "`mask` must have the same shape as `Q` and `U`"
148        
149    def get_first_der(self):
150        raise NotImplementedError("First derivatives are not implemented for arbitrary SO(3) DataFields, please use `SO3Healpix` or `SO3Patch`.")
151    
152    def get_second_der(self):
153        raise NotImplementedError("Second derivatives are not implemented for arbitrary SO(3) DataFields, please use `SO3Healpix` or `SO3Patch`.")
154        
155    def _V0(self, us, dus, verbose=True):
156        """Compute the first Minkowski Functional, $V_0$, normalized by the volume of the space. Internal interface for `pynkowski.V0`.
157
158        Parameters
159        ----------
160        us : np.array
161            The thresholds where $V_0$ is computed.
162
163        dus : np.array
164            The bin sizes of `us`. Ignored.
165
166        verbose : bool, optional
167            If True (default), progress bars are shown for the computations on data.
168        
169        Returns
170        -------
171        V0 : np.array()
172            The values of the first Minkowski Functional at the given thresholds.
173        
174        """
175        stat = np.zeros_like(us)
176        modulus = self.field.modulus()
177
178        for ii in tqdm(np.arange(us.shape[0]), disable=not verbose):
179            lenghts = np.zeros_like(modulus)
180            lenghts[us[ii]<=-modulus] = np.pi
181            mask = (us[ii]>-modulus) & (us[ii]<modulus)
182            lenghts[mask] = np.arccos(us[ii]/modulus[mask])
183            stat[ii] = np.mean(lenghts[self.mask])/np.pi
184        return stat
185    
186    def _V1(self, us, dus, verbose=True):
187        """Compute the second Minkowski Functional, $V_1$, normalized by the volume of the space. Internal interface for `pynkowski.V1`.
188
189        Parameters
190        ----------
191        us : np.array
192            The thresholds where $V_1$ is computed.
193            
194        dus : np.array
195            The bin sizes of `us`. Ignored.
196            
197        verbose : bool, optional
198            If True (default), progress bars are shown for the computations on data.
199        
200        Returns
201        -------
202        V1 : np.array()
203            The values of the second Minkowski Functional at the given thresholds.
204        
205        """
206        if self.first_der is None:
207            self.get_first_der()
208        modulus = self.field.modulus()
209        pol_angle = self.field.pol_angle()
210        # theta = get_theta(self.nside)
211        stat = np.zeros_like(us)
212        
213        for ii in tqdm(np.arange(us.shape[0]), disable=not verbose):
214            integrand_1 = np.zeros_like(modulus)
215            integrand_2 = np.zeros_like(modulus)
216            # thetamask = ~(np.isclose(np.cos(theta),0, atol=1.e-2))
217            mask = (us[ii]>-modulus) & (us[ii]<modulus)   # & thetamask
218            if np.sum(mask) == 0:
219                continue
220            first_der_t = QUarray(*zip(*self.first_der))[:,:,mask]
221
222            psi_1 = pol_angle[mask] + np.arccos(us[ii]/modulus[mask])/2.
223            psi_2 = pol_angle[mask] - np.arccos(us[ii]/modulus[mask])/2.
224            
225            integrand_1[mask] = np.sqrt(np.sum(first_der_t(psi_1)**2., axis=0))
226            integrand_2[mask] = np.sqrt(np.sum(first_der_t(psi_2)**2., axis=0))
227            
228            total_integrand = (integrand_1 + integrand_2)   # /np.mean(thetamask)
229            stat[ii] = np.mean(total_integrand[self.mask])/np.pi
230
231        return _MF_prefactor(self.dim, 1) * stat
232    
233    
234    def _V2(self, us, dus, verbose=True):
235        """Compute the second Minkowski Functional, $V_2$, normalized by the volume of the space. Internal interface for `pynkowski.V2`.
236
237        Parameters
238        ----------
239        us : np.array
240            The thresholds where $V_2$ is computed.
241            
242        dus : np.array
243            The bin sizes of `us`. Ignored.
244            
245        verbose : bool, optional
246            If True (default), progress bars are shown for the computations on data.
247        
248        Returns
249        -------
250        V2 : np.array()
251            The values of the second Minkowski Functional at the given thresholds.
252        
253        """
254        if self.first_der is None:
255            self.get_first_der()
256        if self.second_der is None:
257            self.get_second_der()
258            
259        modulus = self.field.modulus()
260        pol_angle = self.field.pol_angle()
261        # theta = get_theta(self.nside)
262        stat = np.zeros_like(us)
263
264        for ii in tqdm(np.arange(us.shape[0]), disable=not verbose):
265            integrand_1 = np.zeros_like(modulus)
266            integrand_2 = np.zeros_like(modulus)
267            
268            # thetamask = ~(np.isclose(np.cos(theta),0, atol=1.e-2))
269            mask = (us[ii]>-modulus) & (us[ii]<modulus)  #& thetamask
270            if np.sum(mask) == 0:
271                continue
272            # pixs = np.arange(12*self.nside**2)[mask]
273            first_der_t = QUarray(*zip(*self.first_der))[:,:,mask]       # This just "transposes" the array without losing QUarray functionality, so that the first index is the Q/U component, the second is the derivative direction, and the third is the pixel
274            second_der_t = QUarray(*zip(*self.second_der))[:,:,mask]
275            
276            psi_1 = pol_angle[mask] + np.arccos(us[ii]/modulus[mask])/2.
277            psi_2 = pol_angle[mask] - np.arccos(us[ii]/modulus[mask])/2.
278
279            integrand_1[mask] = self._H(first_der_t, second_der_t, psi_1)
280            integrand_2[mask] = self._H(first_der_t, second_der_t, psi_2)
281        
282            total_integrand = (integrand_1 + integrand_2)  #/np.mean(thetamask)
283        
284            stat[ii] = np.mean(total_integrand[self.mask])/np.pi
285
286        return _MF_prefactor(self.dim, 2) * stat
287    
288    def _V3(self, us, dus, verbose=True):
289        """Compute the second Minkowski Functional, $V_3$, normalized by the volume of the space. Internal interface for `pynkowski.V3`.
290
291        Parameters
292        ----------
293        us : np.array
294            The thresholds where $V_3$ is computed.
295            
296        dus : np.array
297            The bin sizes of `us`. Ignored.
298            
299        verbose : bool, optional
300            If True (default), progress bars are shown for the computations on data.
301        
302        Returns
303        -------
304        V3 : np.array()
305            The values of the second Minkowski Functional at the given thresholds.
306        
307        """
308        if self.first_der is None:
309            self.get_first_der()
310        if self.second_der is None:
311            self.get_second_der()
312            
313        modulus = self.field.modulus()
314        pol_angle = self.field.pol_angle()
315        # theta = get_theta(self.nside)
316        stat = np.zeros_like(us)
317
318        for ii in tqdm(np.arange(us.shape[0]), disable=not verbose):
319            integrand_1 = np.zeros_like(modulus)
320            integrand_2 = np.zeros_like(modulus)
321            
322            # thetamask = ~(np.isclose(np.cos(theta),0, atol=1.e-2))
323            mask = (us[ii]>-modulus) & (us[ii]<modulus)  #& thetamask
324            if np.sum(mask) == 0:
325                continue
326            # pixs = np.arange(12*self.nside**2)[mask]
327            first_der_t = QUarray(*zip(*self.first_der))[:,:,mask]
328            second_der_t = QUarray(*zip(*self.second_der))[:,:,mask]
329            
330            psi_1 = pol_angle[mask] + np.arccos(us[ii]/modulus[mask])/2.
331            psi_2 = pol_angle[mask] - np.arccos(us[ii]/modulus[mask])/2.
332
333            integrand_1[mask] = self._K(first_der_t, second_der_t, psi_1)
334            integrand_2[mask] = self._K(first_der_t, second_der_t, psi_2)
335        
336            total_integrand = (integrand_1 + integrand_2)  #/np.mean(thetamask)
337        
338            stat[ii] = np.mean(total_integrand[self.mask])/np.pi
339
340        return _MF_prefactor(self.dim, 3) * stat
341
342    @staticmethod
343    def _H(first_der_t, second_der_t, psi):
344        """Compute the mean curvature (times 2) of the field at a given angle.
345        """
346        first_der_t = first_der_t(psi)
347        second_der_t = second_der_t(psi)
348        hess = np.array([[second_der_t[0], second_der_t[3], second_der_t[4]],
349                        [second_der_t[3], second_der_t[1], second_der_t[5]],
350                        [second_der_t[4], second_der_t[5], second_der_t[2]]])
351        norm_grad = first_der_t / np.sqrt(np.sum(first_der_t**2., axis=0))
352        return np.einsum('j...,jk...,k... -> ...', norm_grad, hess, norm_grad) - np.trace(hess, axis1=0, axis2=1)
353    
354    @staticmethod
355    def _K(first_der_t, second_der_t, psi):
356        """Compute the Gaussian curvature of the field at a given angle.
357        """
358        first_der_t = first_der_t(psi)
359        mod_grad = np.sqrt(np.sum(first_der_t**2., axis=0))
360        second_der_t = second_der_t(psi) / mod_grad
361        norm_grad = first_der_t / mod_grad
362        extended_hessian_det = np.linalg.det(np.array([[second_der_t[0], second_der_t[3], second_der_t[4], norm_grad[0]],
363                                                        [second_der_t[3], second_der_t[1], second_der_t[5], norm_grad[1]],
364                                                        [second_der_t[4], second_der_t[5], second_der_t[2], norm_grad[2]],
365                                                        [norm_grad[0], norm_grad[1], norm_grad[2], np.zeros_like(norm_grad[0])]]).T).T
366        return -extended_hessian_det*mod_grad
367
368
369class SO3Healpix(SO3DataField):
370    """Class for spherical spin fields in the SO(3) formalism, with Q and U in HEALPix format.
371
372    Parameters
373    ----------
374    Q : np.array
375        Values of the Q field in HEALPix format in RING scheme.
376    
377    U : np.array
378        Values of the U field in HEALPix format in RING scheme. Must have the same shape as `Q`.
379        
380    normalise : bool, optional
381        If `True`, the map is normalised to have unit variance.
382        Default: `True`.
383
384    mask : np.array or None, optional
385        Sky mask where the field is considered. It is a bool array of the same shape that `Q` and `U`.
386        Default: all data is included.
387        
388    Attributes
389    ----------
390    field : np.array
391        Array of shape `(2, Q.shape)` containing the Q and U values. It can be called as `field(psi)` to obtain the value of $f$ for the given `psi`.
392        
393    nside : int
394        Parameter `nside` of the Q and U maps. The number of pixels is `12*nside**2`.
395    
396    dim : int
397        Dimension of the space where the field is defined. In this case, the space is SO(3) and this is 3.
398        
399    name : str
400        Name of the field. In this case, it is `"SO(3) HEALPix map"`.
401    
402    first_der : list or None
403        First **covariant** derivatives of the field in an orthonormal basis of the space. It a list of size `3`, each entry has the same structure as `field`.
404    
405    second_der : list or None
406        Second **covariant** derivatives of the field in an orthonormal basis of the space. It a list of size `6`, each entry has the same structure as `field`.
407        The order of the derivatives is diagonal first, i.e., `11`, `22`, `33`, `12`, `13`, `23`.
408        
409    mask : np.array
410        Sky mask where the field is considered. It is a bool array of the same shape that `Q` and `U`.
411
412    """   
413    def __init__(self, Q, U, normalise=True, mask=None):
414        super().__init__(Q, U, name="SO(3) HEALPix map", mask=mask)
415        self.nside = hp.get_nside(self.field[0])
416        if hp.get_nside(self.mask) != self.nside:
417            raise ValueError('The map and the mask have different nside')
418        
419        if normalise:
420            σ2 = self.get_variance()
421            self.field /= np.sqrt(σ2)
422           
423    def get_variance(self):
424        """Compute the variance of the SO(3) Healpix map within the sky mask. 
425
426        Returns
427        -------
428        var : float
429            The variance of the map within the mask.
430
431        """    
432        return (np.mean(self.field[:,self.mask]**2.))
433    
434    def get_first_der(self, lmax=None):
435        """Compute the covariant first derivatives of the SO(3) Healpix map. 
436        It stores:
437        
438        - first covariant derivative wrt e₁ in self.first_der[0]
439        - first covariant derivative wrt e₂ in self.first_der[1]
440        - first covariant derivative wrt e₃ in self.first_der[2]
441        
442        Parameters
443        ----------
444        lmax : int or None, optional
445            Maximum multipole used in the computation of the derivatives.
446            Default: 3*nside - 1
447
448        """    
449        Q_grad = healpix_derivatives(self.field[0], gradient=True, lmax=lmax)  # order θ,ϕ
450        U_grad = healpix_derivatives(self.field[1], gradient=True, lmax=lmax)
451        theta = get_theta(self.nside)
452
453        self._fphi = QUarray(Q_grad[1], U_grad[1])
454        ftheta = QUarray(Q_grad[0], U_grad[0])
455        fpsi = self.field.derpsi()
456
457        factor_p = (np.sqrt(1.-np.sin(theta)) + np.sqrt(1.+np.sin(theta))) / np.sqrt(8.)
458        factor_m = (np.sqrt(1.-np.sin(theta)) - np.sqrt(1.+np.sin(theta))) / np.sqrt(8.)
459
460        self.first_der = [ factor_p * self._fphi + factor_m * fpsi / np.cos(theta),         # ∂e₁
461                           ftheta / np.sqrt(2.),                                            # ∂e₂
462                           factor_m * self._fphi + factor_p * fpsi / np.cos(theta) ]        # ∂e₃
463
464        
465        # This is the old way, I leave it here for reference and debugging if needed
466        # self.grad_theta = Polmap(Q_grad[0]/np.sqrt(2.), U_grad[0]/np.sqrt(2.), normalise=False)
467        # self.grad_phi = Polmap(((np.sqrt(1.-np.sin(theta)) + np.sqrt(1.+np.sin(theta)) ) / np.sqrt(8.) ) * Q_grad[1] -
468        #                           ((np.sqrt(1.-np.sin(theta)) - np.sqrt(1.+np.sin(theta)) ) / (np.sqrt(8.)*np.cos(theta)) ) * 2.*self.U, 
469        #                        ((np.sqrt(1.-np.sin(theta)) + np.sqrt(1.+np.sin(theta)) ) / np.sqrt(8.) ) * U_grad[1] +
470        #                           ((np.sqrt(1.-np.sin(theta)) - np.sqrt(1.+np.sin(theta)) ) / (np.sqrt(8.)*np.cos(theta)) ) * 2.*self.Q, normalise=False)
471        
472        # self.grad_psi = Polmap(((np.sqrt(1.-np.sin(theta)) - np.sqrt(1.+np.sin(theta)) ) / np.sqrt(8.) ) * Q_grad[1] -
473        #                           ((np.sqrt(1.-np.sin(theta)) + np.sqrt(1.+np.sin(theta)) ) / (np.sqrt(8.)*np.cos(theta)) ) * 2.*self.U, 
474        #                        ((np.sqrt(1.-np.sin(theta)) - np.sqrt(1.+np.sin(theta)) ) / np.sqrt(8.) ) * U_grad[1] +
475        #                           ((np.sqrt(1.-np.sin(theta)) + np.sqrt(1.+np.sin(theta)) ) / (np.sqrt(8.)*np.cos(theta)) ) * 2.*self.Q, normalise=False)
476        
477        # self.der_phi = Polmap(np.cos(theta) * Q_grad[1], np.cos(theta) * U_grad[1], normalise=False)
478
479
480        # self.first_der = np.array(healpix_derivatives(self.field, gradient=True, lmax=lmax))  # order: θ, ϕ
481            
482    def get_second_der(self, lmax=None):
483        """Compute the covariant second derivatives of the input Healpix scalar map. 
484        It stores:
485        
486        - second covariant derivative wrt e₁e₁ in self.second_der[0]
487        - second covariant derivative wrt e₂e₂ in self.second_der[1]
488        - second covariant derivative wrt e₃e₃ in self.second_der[2]
489        - second covariant derivative wrt e₁e₂ in self.second_der[3]
490        - second covariant derivative wrt e₁e₃ in self.second_der[4]
491        - second covariant derivative wrt e₂e₃ in self.second_der[5]
492
493        Parameters
494        ----------
495        lmax : int or None, optional
496            Maximum multipole used in the computation of the derivatives.
497            Default: 3*nside - 1
498
499        """    
500        if self.first_der is None:
501            self.get_first_der()
502        theta = get_theta(self.nside)
503
504        # self._fphi  is already computed in get_first_der
505        self._fphi *= np.cos(theta)
506        ftheta = self.first_der[1] * np.sqrt(2.)
507        fpsi = self.field.derpsi()
508
509        Q_der_der = healpix_second_derivatives(ftheta[0], self._fphi[0], lmax=lmax)  # order θθ, ϕϕ, θϕ
510        U_der_der = healpix_second_derivatives(ftheta[1], self._fphi[1], lmax=lmax)
511        fthetatheta = QUarray(Q_der_der[0], U_der_der[0])
512        fphiphi = QUarray(Q_der_der[1], U_der_der[1])
513        fthetaphi = QUarray(Q_der_der[2], U_der_der[2])
514        del Q_der_der, U_der_der
515
516        factor_P = (np.sqrt(1.-np.sin(theta)) + np.sqrt(1.+np.sin(theta))) 
517        factor_M = (np.sqrt(1.-np.sin(theta)) - np.sqrt(1.+np.sin(theta)))
518
519        # order 11, 22, 33, 12, 13, 23
520        self.second_der = [ ((1.+np.cos(theta))*fphiphi + (1.-np.cos(theta))*fpsi.derpsi() - 2.*np.sin(theta)*self._fphi.derpsi() - np.sin(2.*theta)*ftheta/2.) / (4.*np.cos(theta)**2.),
521                            fthetatheta / 2.,
522                            ((1.-np.cos(theta))*fphiphi + (1.+np.cos(theta))*fpsi.derpsi() - 2.*np.sin(theta)*self._fphi.derpsi() - np.sin(2.*theta)*ftheta/2.) / (4.*np.cos(theta)**2.),
523                            (2.*np.cos(theta)*(factor_M*ftheta.derpsi() + factor_P*fthetaphi) + (factor_P*np.sin(theta) - factor_M)*self._fphi + (factor_M*np.sin(theta) - factor_P)*fpsi ) / (8.*np.cos(theta)**2.),
524                            (-np.sin(theta)*fphiphi - np.sin(theta)*fpsi.derpsi() + 2.*self._fphi.derpsi() + np.cos(theta)*ftheta ) / (4.*np.cos(theta)**2.),
525                            (2.*np.cos(theta)*(factor_M*fthetaphi + factor_P*ftheta.derpsi()) + (factor_M*np.sin(theta) - factor_P)*self._fphi + (factor_P*np.sin(theta) - factor_M)*fpsi ) / (8.*np.cos(theta)**2.) ]
526
527        
528        # This is the old way, I leave it here for reference and debugging if needed
529        # self.der_theta_theta = Polmap(Q_der_der[0]/2., U_der_der[0]/2., normalise=False)
530        
531        
532        # self.der_phi_phi = Polmap( ((1.+np.cos(theta)) * Q_der_der[1]  - (1-np.cos(theta)) * 4. * self.Q  + 4.*np.sin(theta)*U_der[1] - np.sin(2.*theta)*Q_der[0]/2. ) / (4.*np.cos(theta)**2.),
533        #                            ((1.+np.cos(theta)) * U_der_der[1]  - (1-np.cos(theta)) * 4. * self.U  - 4.*np.sin(theta)*Q_der[1] - np.sin(2.*theta)*U_der[0]/2. ) / (4.*np.cos(theta)**2.), normalise=False)
534            
535        # self.der_psi_psi = Polmap( ((1.-np.cos(theta)) * Q_der_der[1]  - (1+np.cos(theta)) * 4. * self.Q  + 4.*np.sin(theta)*U_der[1] - np.sin(2.*theta)*Q_der[0]/2. ) / (4.*np.cos(theta)**2.),
536        #                            ((1.-np.cos(theta)) * U_der_der[1]  - (1+np.cos(theta)) * 4. * self.U  - 4.*np.sin(theta)*Q_der[1] - np.sin(2.*theta)*U_der[0]/2. ) / (4.*np.cos(theta)**2.), normalise=False)
537            
538        # P = np.sqrt(1.-np.sin(theta)) + np.sqrt(1.+np.sin(theta))
539        # M = np.sqrt(1.-np.sin(theta)) - np.sqrt(1.+np.sin(theta))
540        
541        # self.der_theta_phi = Polmap( (2.*np.cos(theta)* (-2.*M*U_der[0] + P*Q_der_der[2]) + (P*np.sin(theta)-M)*Q_der[1] - 2.*(M*np.sin(theta)-P)*self.U ) / (8.*np.cos(theta)**2.),
542        #                              (2.*np.cos(theta)* ( 2.*M*Q_der[0] + P*U_der_der[2]) + (P*np.sin(theta)-M)*U_der[1] + 2.*(M*np.sin(theta)-P)*self.Q ) / (8.*np.cos(theta)**2.), normalise=False)
543            
544        # self.der_theta_psi = Polmap( (2.*np.cos(theta)* (-2.*P*U_der[0] + M*Q_der_der[2]) + (M*np.sin(theta)-P)*Q_der[1] - 2.*(P*np.sin(theta)-M)*self.U ) / (8.*np.cos(theta)**2.),
545        #                              (2.*np.cos(theta)* ( 2.*P*Q_der[0] + M*U_der_der[2]) + (M*np.sin(theta)-P)*U_der[1] + 2.*(P*np.sin(theta)-M)*self.Q ) / (8.*np.cos(theta)**2.), normalise=False)
546        # # self.der_theta_psi = Polmap( (2.*np.cos(theta)* (-2.*P*U_der[0] + M*Q_der_der[2]) + ( 3.*M*np.sin(theta) - P*(1.+np.cos(theta)) )*Q_der[1] - 2.*(3.*P*np.sin(theta) - M*(1.-np.cos(theta)))*self.U ) / (8.*np.cos(theta)**2.),
547        # #                              (2.*np.cos(theta)* ( 2.*P*Q_der[0] + M*U_der_der[2]) + ( 3.*M*np.sin(theta) - P*(1.+np.cos(theta)) )*U_der[1] + 2.*(3.*P*np.sin(theta) - M*(1.-np.cos(theta)))*self.Q ) / (8.*np.cos(theta)**2.), normalise=False)
548            
549        # self.der_phi_psi = Polmap( (-np.sin(theta) * Q_der_der[1]  + np.sin(theta) * 4. * self.Q  - 4.*U_der[1] + np.cos(theta)*Q_der[0] ) / (4.*np.cos(theta)**2.),
550        #                            (-np.sin(theta) * U_der_der[1]  + np.sin(theta) * 4. * self.U  + 4.*Q_der[1] + np.cos(theta)*U_der[0] ) / (4.*np.cos(theta)**2.), normalise=False)
551
552
553class SO3Patch(SO3DataField):
554    """Class for spin fields in the SO(3) formalism, with Q and U as 2D `np.array` patches. This assumes a flat sky approximation on the patch.
555
556    Parameters
557    ----------
558    Q : np.array
559        Values of the Q field as a 2D array. The x, y dimensions correspond to the ϕ, θ coordinates.
560    
561    U : np.array
562        Values of the U field as a 2D array. Must have the same shape as `Q`.
563        
564    normalise : bool, optional
565        If `True`, the map is normalised to have unit variance.
566        Default: `True`.
567
568    mask : np.array or None, optional
569        Mask where the field is considered. It is a bool array of the same shape that `Q` and `U`.
570        Default: all data is included.
571
572    spacing : float or np.array, optional
573        Spacing between pixels (centres) in each dimension, in radians. If a float, the spacing is the same in all dimensions. 
574        If an array, it must have the same length as the number of dimensions of the field (2).
575        
576    Attributes
577    ----------
578    field : np.array
579        Array of shape `(2, Q.shape)` containing the Q and U values. It can be called as `field(psi)` to obtain the value of $f$ for the given `psi`.
580            
581    dim : int
582        Dimension of the space where the field is defined. In this case, the space is SO(3) and this is 3.
583        
584    name : str
585        Name of the field. In this case, it is `"SO(3) patch"`.
586    
587    first_der : list or None
588        First **covariant** derivatives of the field in an orthonormal basis of the space. It a list of size `3`, each entry has the same structure as `field`.
589    
590    second_der : list or None
591        Second **covariant** derivatives of the field in an orthonormal basis of the space. It a list of size `6`, each entry has the same structure as `field`.
592        The order of the derivatives is diagonal first, i.e., `11`, `22`, `33`, `12`, `13`, `23`.
593        
594    mask : np.array
595        Mask where the field is considered. It is a bool array of the same shape that `Q` and `U`.
596
597    spacing : float or np.array
598        Spacing between pixels (centres) in each dimension. If a float, the spacing is the same in all dimensions.
599        If an array, it must have the same length as the number of dimensions of the field (2).
600
601    """   
602    def __init__(self, Q, U, normalise=True, mask=None, spacing=np.deg2rad(1./60.)):
603        super().__init__(Q, U, name="SO(3) HEALPix map", mask=mask)
604        if self.mask.shape != Q.shape:
605            raise ValueError('The map and the mask have different shapes')
606        
607        if normalise:
608            σ2 = self.get_variance()
609            self.field /= np.sqrt(σ2)
610        self.spacing = spacing
611           
612    def get_variance(self):
613        """Compute the variance of the SO(3) Healpix map within the sky mask. 
614
615        Returns
616        -------
617        var : float
618            The variance of the map within the mask.
619
620        """    
621        return (np.mean(self.field[:,self.mask]**2.))
622    
623    def get_first_der(self):
624        """Compute the covariant first derivatives of the SO(3) Healpix map. Flat sky is assumed.
625        It stores:
626        
627        - first covariant derivative wrt ϕ in self.first_der[0]
628        - first covariant derivative wrt θ in self.first_der[1]
629        - first covariant derivative wrt ψ in self.first_der[2]
630        
631        """    
632        Q_grad = np.gradient(self.field[0], self.spacing, edge_order=2)[::-1]  # order θ,ϕ
633        U_grad = np.gradient(self.field[1], self.spacing, edge_order=2)[::-1]
634
635        fphi = QUarray(Q_grad[1], U_grad[1])
636        ftheta = QUarray(Q_grad[0], U_grad[0])
637        fpsi = self.field.derpsi()
638
639        self.first_der = [ fphi / np.sqrt(2.),                  # ∂e₁  (ϕ)
640                           ftheta / np.sqrt(2.),                # ∂e₂  (θ)
641                           fpsi / np.sqrt(2.)]                  # ∂e₃  (ψ)
642
643        
644            
645    def get_second_der(self):
646        """Compute the covariant second derivatives of the input Healpix scalar map. 
647        It stores:
648        
649        - second covariant derivative wrt ϕϕ in self.second_der[0]
650        - second covariant derivative wrt θθ in self.second_der[1]
651        - second covariant derivative wrt ψψ in self.second_der[2]
652        - second covariant derivative wrt ϕθ in self.second_der[3]
653        - second covariant derivative wrt ϕψ in self.second_der[4]
654        - second covariant derivative wrt θψ in self.second_der[5]
655
656        """    
657        if self.first_der is None:
658            self.get_first_der()
659
660        fphi = self.first_der[0] * np.sqrt(2.)
661        ftheta = self.first_der[1] * np.sqrt(2.)
662        fpsi = self.field.derpsi()
663
664        Q_phi_der = np.gradient(fphi[0], self.spacing, edge_order=2)[::-1]  # order θ,ϕ
665        U_phi_der = np.gradient(fphi[1], self.spacing, edge_order=2)[::-1]  # order θ,ϕ
666        Q_the_der = np.gradient(ftheta[0], self.spacing, edge_order=2)[::-1]
667        U_the_der = np.gradient(ftheta[1], self.spacing, edge_order=2)[::-1]
668
669        fthetatheta = QUarray(Q_the_der[0], U_the_der[0])
670        fthetaphi = QUarray(Q_the_der[1], U_the_der[1])
671        fphiphi = QUarray(Q_phi_der[1], U_phi_der[1])
672        del Q_the_der, U_the_der, Q_phi_der, U_phi_der
673
674        # order 11, 22, 33, 12, 13, 23
675        self.second_der = [ fphiphi / 2.,
676                            fthetatheta / 2.,
677                            fpsi.derpsi() / 2.,
678                            fthetaphi / 2.,
679                            fphi.derpsi() / 2.,
680                            ftheta.derpsi() / 2. ]
681
682
683__all__ = ["SO3Healpix", "SO3Patch", "QUarray", "SO3DataField"]
684
685__docformat__ = "numpy"
class SO3Healpix(SO3DataField):
370class SO3Healpix(SO3DataField):
371    """Class for spherical spin fields in the SO(3) formalism, with Q and U in HEALPix format.
372
373    Parameters
374    ----------
375    Q : np.array
376        Values of the Q field in HEALPix format in RING scheme.
377    
378    U : np.array
379        Values of the U field in HEALPix format in RING scheme. Must have the same shape as `Q`.
380        
381    normalise : bool, optional
382        If `True`, the map is normalised to have unit variance.
383        Default: `True`.
384
385    mask : np.array or None, optional
386        Sky mask where the field is considered. It is a bool array of the same shape that `Q` and `U`.
387        Default: all data is included.
388        
389    Attributes
390    ----------
391    field : np.array
392        Array of shape `(2, Q.shape)` containing the Q and U values. It can be called as `field(psi)` to obtain the value of $f$ for the given `psi`.
393        
394    nside : int
395        Parameter `nside` of the Q and U maps. The number of pixels is `12*nside**2`.
396    
397    dim : int
398        Dimension of the space where the field is defined. In this case, the space is SO(3) and this is 3.
399        
400    name : str
401        Name of the field. In this case, it is `"SO(3) HEALPix map"`.
402    
403    first_der : list or None
404        First **covariant** derivatives of the field in an orthonormal basis of the space. It a list of size `3`, each entry has the same structure as `field`.
405    
406    second_der : list or None
407        Second **covariant** derivatives of the field in an orthonormal basis of the space. It a list of size `6`, each entry has the same structure as `field`.
408        The order of the derivatives is diagonal first, i.e., `11`, `22`, `33`, `12`, `13`, `23`.
409        
410    mask : np.array
411        Sky mask where the field is considered. It is a bool array of the same shape that `Q` and `U`.
412
413    """   
414    def __init__(self, Q, U, normalise=True, mask=None):
415        super().__init__(Q, U, name="SO(3) HEALPix map", mask=mask)
416        self.nside = hp.get_nside(self.field[0])
417        if hp.get_nside(self.mask) != self.nside:
418            raise ValueError('The map and the mask have different nside')
419        
420        if normalise:
421            σ2 = self.get_variance()
422            self.field /= np.sqrt(σ2)
423           
424    def get_variance(self):
425        """Compute the variance of the SO(3) Healpix map within the sky mask. 
426
427        Returns
428        -------
429        var : float
430            The variance of the map within the mask.
431
432        """    
433        return (np.mean(self.field[:,self.mask]**2.))
434    
435    def get_first_der(self, lmax=None):
436        """Compute the covariant first derivatives of the SO(3) Healpix map. 
437        It stores:
438        
439        - first covariant derivative wrt e₁ in self.first_der[0]
440        - first covariant derivative wrt e₂ in self.first_der[1]
441        - first covariant derivative wrt e₃ in self.first_der[2]
442        
443        Parameters
444        ----------
445        lmax : int or None, optional
446            Maximum multipole used in the computation of the derivatives.
447            Default: 3*nside - 1
448
449        """    
450        Q_grad = healpix_derivatives(self.field[0], gradient=True, lmax=lmax)  # order θ,ϕ
451        U_grad = healpix_derivatives(self.field[1], gradient=True, lmax=lmax)
452        theta = get_theta(self.nside)
453
454        self._fphi = QUarray(Q_grad[1], U_grad[1])
455        ftheta = QUarray(Q_grad[0], U_grad[0])
456        fpsi = self.field.derpsi()
457
458        factor_p = (np.sqrt(1.-np.sin(theta)) + np.sqrt(1.+np.sin(theta))) / np.sqrt(8.)
459        factor_m = (np.sqrt(1.-np.sin(theta)) - np.sqrt(1.+np.sin(theta))) / np.sqrt(8.)
460
461        self.first_der = [ factor_p * self._fphi + factor_m * fpsi / np.cos(theta),         # ∂e₁
462                           ftheta / np.sqrt(2.),                                            # ∂e₂
463                           factor_m * self._fphi + factor_p * fpsi / np.cos(theta) ]        # ∂e₃
464
465        
466        # This is the old way, I leave it here for reference and debugging if needed
467        # self.grad_theta = Polmap(Q_grad[0]/np.sqrt(2.), U_grad[0]/np.sqrt(2.), normalise=False)
468        # self.grad_phi = Polmap(((np.sqrt(1.-np.sin(theta)) + np.sqrt(1.+np.sin(theta)) ) / np.sqrt(8.) ) * Q_grad[1] -
469        #                           ((np.sqrt(1.-np.sin(theta)) - np.sqrt(1.+np.sin(theta)) ) / (np.sqrt(8.)*np.cos(theta)) ) * 2.*self.U, 
470        #                        ((np.sqrt(1.-np.sin(theta)) + np.sqrt(1.+np.sin(theta)) ) / np.sqrt(8.) ) * U_grad[1] +
471        #                           ((np.sqrt(1.-np.sin(theta)) - np.sqrt(1.+np.sin(theta)) ) / (np.sqrt(8.)*np.cos(theta)) ) * 2.*self.Q, normalise=False)
472        
473        # self.grad_psi = Polmap(((np.sqrt(1.-np.sin(theta)) - np.sqrt(1.+np.sin(theta)) ) / np.sqrt(8.) ) * Q_grad[1] -
474        #                           ((np.sqrt(1.-np.sin(theta)) + np.sqrt(1.+np.sin(theta)) ) / (np.sqrt(8.)*np.cos(theta)) ) * 2.*self.U, 
475        #                        ((np.sqrt(1.-np.sin(theta)) - np.sqrt(1.+np.sin(theta)) ) / np.sqrt(8.) ) * U_grad[1] +
476        #                           ((np.sqrt(1.-np.sin(theta)) + np.sqrt(1.+np.sin(theta)) ) / (np.sqrt(8.)*np.cos(theta)) ) * 2.*self.Q, normalise=False)
477        
478        # self.der_phi = Polmap(np.cos(theta) * Q_grad[1], np.cos(theta) * U_grad[1], normalise=False)
479
480
481        # self.first_der = np.array(healpix_derivatives(self.field, gradient=True, lmax=lmax))  # order: θ, ϕ
482            
483    def get_second_der(self, lmax=None):
484        """Compute the covariant second derivatives of the input Healpix scalar map. 
485        It stores:
486        
487        - second covariant derivative wrt e₁e₁ in self.second_der[0]
488        - second covariant derivative wrt e₂e₂ in self.second_der[1]
489        - second covariant derivative wrt e₃e₃ in self.second_der[2]
490        - second covariant derivative wrt e₁e₂ in self.second_der[3]
491        - second covariant derivative wrt e₁e₃ in self.second_der[4]
492        - second covariant derivative wrt e₂e₃ in self.second_der[5]
493
494        Parameters
495        ----------
496        lmax : int or None, optional
497            Maximum multipole used in the computation of the derivatives.
498            Default: 3*nside - 1
499
500        """    
501        if self.first_der is None:
502            self.get_first_der()
503        theta = get_theta(self.nside)
504
505        # self._fphi  is already computed in get_first_der
506        self._fphi *= np.cos(theta)
507        ftheta = self.first_der[1] * np.sqrt(2.)
508        fpsi = self.field.derpsi()
509
510        Q_der_der = healpix_second_derivatives(ftheta[0], self._fphi[0], lmax=lmax)  # order θθ, ϕϕ, θϕ
511        U_der_der = healpix_second_derivatives(ftheta[1], self._fphi[1], lmax=lmax)
512        fthetatheta = QUarray(Q_der_der[0], U_der_der[0])
513        fphiphi = QUarray(Q_der_der[1], U_der_der[1])
514        fthetaphi = QUarray(Q_der_der[2], U_der_der[2])
515        del Q_der_der, U_der_der
516
517        factor_P = (np.sqrt(1.-np.sin(theta)) + np.sqrt(1.+np.sin(theta))) 
518        factor_M = (np.sqrt(1.-np.sin(theta)) - np.sqrt(1.+np.sin(theta)))
519
520        # order 11, 22, 33, 12, 13, 23
521        self.second_der = [ ((1.+np.cos(theta))*fphiphi + (1.-np.cos(theta))*fpsi.derpsi() - 2.*np.sin(theta)*self._fphi.derpsi() - np.sin(2.*theta)*ftheta/2.) / (4.*np.cos(theta)**2.),
522                            fthetatheta / 2.,
523                            ((1.-np.cos(theta))*fphiphi + (1.+np.cos(theta))*fpsi.derpsi() - 2.*np.sin(theta)*self._fphi.derpsi() - np.sin(2.*theta)*ftheta/2.) / (4.*np.cos(theta)**2.),
524                            (2.*np.cos(theta)*(factor_M*ftheta.derpsi() + factor_P*fthetaphi) + (factor_P*np.sin(theta) - factor_M)*self._fphi + (factor_M*np.sin(theta) - factor_P)*fpsi ) / (8.*np.cos(theta)**2.),
525                            (-np.sin(theta)*fphiphi - np.sin(theta)*fpsi.derpsi() + 2.*self._fphi.derpsi() + np.cos(theta)*ftheta ) / (4.*np.cos(theta)**2.),
526                            (2.*np.cos(theta)*(factor_M*fthetaphi + factor_P*ftheta.derpsi()) + (factor_M*np.sin(theta) - factor_P)*self._fphi + (factor_P*np.sin(theta) - factor_M)*fpsi ) / (8.*np.cos(theta)**2.) ]
527
528        
529        # This is the old way, I leave it here for reference and debugging if needed
530        # self.der_theta_theta = Polmap(Q_der_der[0]/2., U_der_der[0]/2., normalise=False)
531        
532        
533        # self.der_phi_phi = Polmap( ((1.+np.cos(theta)) * Q_der_der[1]  - (1-np.cos(theta)) * 4. * self.Q  + 4.*np.sin(theta)*U_der[1] - np.sin(2.*theta)*Q_der[0]/2. ) / (4.*np.cos(theta)**2.),
534        #                            ((1.+np.cos(theta)) * U_der_der[1]  - (1-np.cos(theta)) * 4. * self.U  - 4.*np.sin(theta)*Q_der[1] - np.sin(2.*theta)*U_der[0]/2. ) / (4.*np.cos(theta)**2.), normalise=False)
535            
536        # self.der_psi_psi = Polmap( ((1.-np.cos(theta)) * Q_der_der[1]  - (1+np.cos(theta)) * 4. * self.Q  + 4.*np.sin(theta)*U_der[1] - np.sin(2.*theta)*Q_der[0]/2. ) / (4.*np.cos(theta)**2.),
537        #                            ((1.-np.cos(theta)) * U_der_der[1]  - (1+np.cos(theta)) * 4. * self.U  - 4.*np.sin(theta)*Q_der[1] - np.sin(2.*theta)*U_der[0]/2. ) / (4.*np.cos(theta)**2.), normalise=False)
538            
539        # P = np.sqrt(1.-np.sin(theta)) + np.sqrt(1.+np.sin(theta))
540        # M = np.sqrt(1.-np.sin(theta)) - np.sqrt(1.+np.sin(theta))
541        
542        # self.der_theta_phi = Polmap( (2.*np.cos(theta)* (-2.*M*U_der[0] + P*Q_der_der[2]) + (P*np.sin(theta)-M)*Q_der[1] - 2.*(M*np.sin(theta)-P)*self.U ) / (8.*np.cos(theta)**2.),
543        #                              (2.*np.cos(theta)* ( 2.*M*Q_der[0] + P*U_der_der[2]) + (P*np.sin(theta)-M)*U_der[1] + 2.*(M*np.sin(theta)-P)*self.Q ) / (8.*np.cos(theta)**2.), normalise=False)
544            
545        # self.der_theta_psi = Polmap( (2.*np.cos(theta)* (-2.*P*U_der[0] + M*Q_der_der[2]) + (M*np.sin(theta)-P)*Q_der[1] - 2.*(P*np.sin(theta)-M)*self.U ) / (8.*np.cos(theta)**2.),
546        #                              (2.*np.cos(theta)* ( 2.*P*Q_der[0] + M*U_der_der[2]) + (M*np.sin(theta)-P)*U_der[1] + 2.*(P*np.sin(theta)-M)*self.Q ) / (8.*np.cos(theta)**2.), normalise=False)
547        # # self.der_theta_psi = Polmap( (2.*np.cos(theta)* (-2.*P*U_der[0] + M*Q_der_der[2]) + ( 3.*M*np.sin(theta) - P*(1.+np.cos(theta)) )*Q_der[1] - 2.*(3.*P*np.sin(theta) - M*(1.-np.cos(theta)))*self.U ) / (8.*np.cos(theta)**2.),
548        # #                              (2.*np.cos(theta)* ( 2.*P*Q_der[0] + M*U_der_der[2]) + ( 3.*M*np.sin(theta) - P*(1.+np.cos(theta)) )*U_der[1] + 2.*(3.*P*np.sin(theta) - M*(1.-np.cos(theta)))*self.Q ) / (8.*np.cos(theta)**2.), normalise=False)
549            
550        # self.der_phi_psi = Polmap( (-np.sin(theta) * Q_der_der[1]  + np.sin(theta) * 4. * self.Q  - 4.*U_der[1] + np.cos(theta)*Q_der[0] ) / (4.*np.cos(theta)**2.),
551        #                            (-np.sin(theta) * U_der_der[1]  + np.sin(theta) * 4. * self.U  + 4.*Q_der[1] + np.cos(theta)*U_der[0] ) / (4.*np.cos(theta)**2.), normalise=False)

Class for spherical spin fields in the SO(3) formalism, with Q and U in HEALPix format.

Parameters
  • Q (np.array): Values of the Q field in HEALPix format in RING scheme.
  • U (np.array): Values of the U field in HEALPix format in RING scheme. Must have the same shape as Q.
  • normalise (bool, optional): If True, the map is normalised to have unit variance. Default: True.
  • mask (np.array or None, optional): Sky mask where the field is considered. It is a bool array of the same shape that Q and U. Default: all data is included.
Attributes
  • field (np.array): Array of shape (2, Q.shape) containing the Q and U values. It can be called as field(psi) to obtain the value of $f$ for the given psi.
  • nside (int): Parameter nside of the Q and U maps. The number of pixels is 12*nside**2.
  • dim (int): Dimension of the space where the field is defined. In this case, the space is SO(3) and this is 3.
  • name (str): Name of the field. In this case, it is "SO(3) HEALPix map".
  • first_der (list or None): First covariant derivatives of the field in an orthonormal basis of the space. It a list of size 3, each entry has the same structure as field.
  • second_der (list or None): Second covariant derivatives of the field in an orthonormal basis of the space. It a list of size 6, each entry has the same structure as field. The order of the derivatives is diagonal first, i.e., 11, 22, 33, 12, 13, 23.
  • mask (np.array): Sky mask where the field is considered. It is a bool array of the same shape that Q and U.
SO3Healpix(Q, U, normalise=True, mask=None)
414    def __init__(self, Q, U, normalise=True, mask=None):
415        super().__init__(Q, U, name="SO(3) HEALPix map", mask=mask)
416        self.nside = hp.get_nside(self.field[0])
417        if hp.get_nside(self.mask) != self.nside:
418            raise ValueError('The map and the mask have different nside')
419        
420        if normalise:
421            σ2 = self.get_variance()
422            self.field /= np.sqrt(σ2)
def get_variance(self):
424    def get_variance(self):
425        """Compute the variance of the SO(3) Healpix map within the sky mask. 
426
427        Returns
428        -------
429        var : float
430            The variance of the map within the mask.
431
432        """    
433        return (np.mean(self.field[:,self.mask]**2.))

Compute the variance of the SO(3) Healpix map within the sky mask.

Returns
  • var (float): The variance of the map within the mask.
def get_first_der(self, lmax=None):
435    def get_first_der(self, lmax=None):
436        """Compute the covariant first derivatives of the SO(3) Healpix map. 
437        It stores:
438        
439        - first covariant derivative wrt e₁ in self.first_der[0]
440        - first covariant derivative wrt e₂ in self.first_der[1]
441        - first covariant derivative wrt e₃ in self.first_der[2]
442        
443        Parameters
444        ----------
445        lmax : int or None, optional
446            Maximum multipole used in the computation of the derivatives.
447            Default: 3*nside - 1
448
449        """    
450        Q_grad = healpix_derivatives(self.field[0], gradient=True, lmax=lmax)  # order θ,ϕ
451        U_grad = healpix_derivatives(self.field[1], gradient=True, lmax=lmax)
452        theta = get_theta(self.nside)
453
454        self._fphi = QUarray(Q_grad[1], U_grad[1])
455        ftheta = QUarray(Q_grad[0], U_grad[0])
456        fpsi = self.field.derpsi()
457
458        factor_p = (np.sqrt(1.-np.sin(theta)) + np.sqrt(1.+np.sin(theta))) / np.sqrt(8.)
459        factor_m = (np.sqrt(1.-np.sin(theta)) - np.sqrt(1.+np.sin(theta))) / np.sqrt(8.)
460
461        self.first_der = [ factor_p * self._fphi + factor_m * fpsi / np.cos(theta),         # ∂e₁
462                           ftheta / np.sqrt(2.),                                            # ∂e₂
463                           factor_m * self._fphi + factor_p * fpsi / np.cos(theta) ]        # ∂e₃
464
465        
466        # This is the old way, I leave it here for reference and debugging if needed
467        # self.grad_theta = Polmap(Q_grad[0]/np.sqrt(2.), U_grad[0]/np.sqrt(2.), normalise=False)
468        # self.grad_phi = Polmap(((np.sqrt(1.-np.sin(theta)) + np.sqrt(1.+np.sin(theta)) ) / np.sqrt(8.) ) * Q_grad[1] -
469        #                           ((np.sqrt(1.-np.sin(theta)) - np.sqrt(1.+np.sin(theta)) ) / (np.sqrt(8.)*np.cos(theta)) ) * 2.*self.U, 
470        #                        ((np.sqrt(1.-np.sin(theta)) + np.sqrt(1.+np.sin(theta)) ) / np.sqrt(8.) ) * U_grad[1] +
471        #                           ((np.sqrt(1.-np.sin(theta)) - np.sqrt(1.+np.sin(theta)) ) / (np.sqrt(8.)*np.cos(theta)) ) * 2.*self.Q, normalise=False)
472        
473        # self.grad_psi = Polmap(((np.sqrt(1.-np.sin(theta)) - np.sqrt(1.+np.sin(theta)) ) / np.sqrt(8.) ) * Q_grad[1] -
474        #                           ((np.sqrt(1.-np.sin(theta)) + np.sqrt(1.+np.sin(theta)) ) / (np.sqrt(8.)*np.cos(theta)) ) * 2.*self.U, 
475        #                        ((np.sqrt(1.-np.sin(theta)) - np.sqrt(1.+np.sin(theta)) ) / np.sqrt(8.) ) * U_grad[1] +
476        #                           ((np.sqrt(1.-np.sin(theta)) + np.sqrt(1.+np.sin(theta)) ) / (np.sqrt(8.)*np.cos(theta)) ) * 2.*self.Q, normalise=False)
477        
478        # self.der_phi = Polmap(np.cos(theta) * Q_grad[1], np.cos(theta) * U_grad[1], normalise=False)
479
480
481        # self.first_der = np.array(healpix_derivatives(self.field, gradient=True, lmax=lmax))  # order: θ, ϕ

Compute the covariant first derivatives of the SO(3) Healpix map. It stores:

  • first covariant derivative wrt e₁ in self.first_der[0]
  • first covariant derivative wrt e₂ in self.first_der[1]
  • first covariant derivative wrt e₃ in self.first_der[2]
Parameters
  • lmax (int or None, optional): Maximum multipole used in the computation of the derivatives. Default: 3*nside - 1
def get_second_der(self, lmax=None):
483    def get_second_der(self, lmax=None):
484        """Compute the covariant second derivatives of the input Healpix scalar map. 
485        It stores:
486        
487        - second covariant derivative wrt e₁e₁ in self.second_der[0]
488        - second covariant derivative wrt e₂e₂ in self.second_der[1]
489        - second covariant derivative wrt e₃e₃ in self.second_der[2]
490        - second covariant derivative wrt e₁e₂ in self.second_der[3]
491        - second covariant derivative wrt e₁e₃ in self.second_der[4]
492        - second covariant derivative wrt e₂e₃ in self.second_der[5]
493
494        Parameters
495        ----------
496        lmax : int or None, optional
497            Maximum multipole used in the computation of the derivatives.
498            Default: 3*nside - 1
499
500        """    
501        if self.first_der is None:
502            self.get_first_der()
503        theta = get_theta(self.nside)
504
505        # self._fphi  is already computed in get_first_der
506        self._fphi *= np.cos(theta)
507        ftheta = self.first_der[1] * np.sqrt(2.)
508        fpsi = self.field.derpsi()
509
510        Q_der_der = healpix_second_derivatives(ftheta[0], self._fphi[0], lmax=lmax)  # order θθ, ϕϕ, θϕ
511        U_der_der = healpix_second_derivatives(ftheta[1], self._fphi[1], lmax=lmax)
512        fthetatheta = QUarray(Q_der_der[0], U_der_der[0])
513        fphiphi = QUarray(Q_der_der[1], U_der_der[1])
514        fthetaphi = QUarray(Q_der_der[2], U_der_der[2])
515        del Q_der_der, U_der_der
516
517        factor_P = (np.sqrt(1.-np.sin(theta)) + np.sqrt(1.+np.sin(theta))) 
518        factor_M = (np.sqrt(1.-np.sin(theta)) - np.sqrt(1.+np.sin(theta)))
519
520        # order 11, 22, 33, 12, 13, 23
521        self.second_der = [ ((1.+np.cos(theta))*fphiphi + (1.-np.cos(theta))*fpsi.derpsi() - 2.*np.sin(theta)*self._fphi.derpsi() - np.sin(2.*theta)*ftheta/2.) / (4.*np.cos(theta)**2.),
522                            fthetatheta / 2.,
523                            ((1.-np.cos(theta))*fphiphi + (1.+np.cos(theta))*fpsi.derpsi() - 2.*np.sin(theta)*self._fphi.derpsi() - np.sin(2.*theta)*ftheta/2.) / (4.*np.cos(theta)**2.),
524                            (2.*np.cos(theta)*(factor_M*ftheta.derpsi() + factor_P*fthetaphi) + (factor_P*np.sin(theta) - factor_M)*self._fphi + (factor_M*np.sin(theta) - factor_P)*fpsi ) / (8.*np.cos(theta)**2.),
525                            (-np.sin(theta)*fphiphi - np.sin(theta)*fpsi.derpsi() + 2.*self._fphi.derpsi() + np.cos(theta)*ftheta ) / (4.*np.cos(theta)**2.),
526                            (2.*np.cos(theta)*(factor_M*fthetaphi + factor_P*ftheta.derpsi()) + (factor_M*np.sin(theta) - factor_P)*self._fphi + (factor_P*np.sin(theta) - factor_M)*fpsi ) / (8.*np.cos(theta)**2.) ]
527
528        
529        # This is the old way, I leave it here for reference and debugging if needed
530        # self.der_theta_theta = Polmap(Q_der_der[0]/2., U_der_der[0]/2., normalise=False)
531        
532        
533        # self.der_phi_phi = Polmap( ((1.+np.cos(theta)) * Q_der_der[1]  - (1-np.cos(theta)) * 4. * self.Q  + 4.*np.sin(theta)*U_der[1] - np.sin(2.*theta)*Q_der[0]/2. ) / (4.*np.cos(theta)**2.),
534        #                            ((1.+np.cos(theta)) * U_der_der[1]  - (1-np.cos(theta)) * 4. * self.U  - 4.*np.sin(theta)*Q_der[1] - np.sin(2.*theta)*U_der[0]/2. ) / (4.*np.cos(theta)**2.), normalise=False)
535            
536        # self.der_psi_psi = Polmap( ((1.-np.cos(theta)) * Q_der_der[1]  - (1+np.cos(theta)) * 4. * self.Q  + 4.*np.sin(theta)*U_der[1] - np.sin(2.*theta)*Q_der[0]/2. ) / (4.*np.cos(theta)**2.),
537        #                            ((1.-np.cos(theta)) * U_der_der[1]  - (1+np.cos(theta)) * 4. * self.U  - 4.*np.sin(theta)*Q_der[1] - np.sin(2.*theta)*U_der[0]/2. ) / (4.*np.cos(theta)**2.), normalise=False)
538            
539        # P = np.sqrt(1.-np.sin(theta)) + np.sqrt(1.+np.sin(theta))
540        # M = np.sqrt(1.-np.sin(theta)) - np.sqrt(1.+np.sin(theta))
541        
542        # self.der_theta_phi = Polmap( (2.*np.cos(theta)* (-2.*M*U_der[0] + P*Q_der_der[2]) + (P*np.sin(theta)-M)*Q_der[1] - 2.*(M*np.sin(theta)-P)*self.U ) / (8.*np.cos(theta)**2.),
543        #                              (2.*np.cos(theta)* ( 2.*M*Q_der[0] + P*U_der_der[2]) + (P*np.sin(theta)-M)*U_der[1] + 2.*(M*np.sin(theta)-P)*self.Q ) / (8.*np.cos(theta)**2.), normalise=False)
544            
545        # self.der_theta_psi = Polmap( (2.*np.cos(theta)* (-2.*P*U_der[0] + M*Q_der_der[2]) + (M*np.sin(theta)-P)*Q_der[1] - 2.*(P*np.sin(theta)-M)*self.U ) / (8.*np.cos(theta)**2.),
546        #                              (2.*np.cos(theta)* ( 2.*P*Q_der[0] + M*U_der_der[2]) + (M*np.sin(theta)-P)*U_der[1] + 2.*(P*np.sin(theta)-M)*self.Q ) / (8.*np.cos(theta)**2.), normalise=False)
547        # # self.der_theta_psi = Polmap( (2.*np.cos(theta)* (-2.*P*U_der[0] + M*Q_der_der[2]) + ( 3.*M*np.sin(theta) - P*(1.+np.cos(theta)) )*Q_der[1] - 2.*(3.*P*np.sin(theta) - M*(1.-np.cos(theta)))*self.U ) / (8.*np.cos(theta)**2.),
548        # #                              (2.*np.cos(theta)* ( 2.*P*Q_der[0] + M*U_der_der[2]) + ( 3.*M*np.sin(theta) - P*(1.+np.cos(theta)) )*U_der[1] + 2.*(3.*P*np.sin(theta) - M*(1.-np.cos(theta)))*self.Q ) / (8.*np.cos(theta)**2.), normalise=False)
549            
550        # self.der_phi_psi = Polmap( (-np.sin(theta) * Q_der_der[1]  + np.sin(theta) * 4. * self.Q  - 4.*U_der[1] + np.cos(theta)*Q_der[0] ) / (4.*np.cos(theta)**2.),
551        #                            (-np.sin(theta) * U_der_der[1]  + np.sin(theta) * 4. * self.U  + 4.*Q_der[1] + np.cos(theta)*U_der[0] ) / (4.*np.cos(theta)**2.), normalise=False)

Compute the covariant second derivatives of the input Healpix scalar map. It stores:

  • second covariant derivative wrt e₁e₁ in self.second_der[0]
  • second covariant derivative wrt e₂e₂ in self.second_der[1]
  • second covariant derivative wrt e₃e₃ in self.second_der[2]
  • second covariant derivative wrt e₁e₂ in self.second_der[3]
  • second covariant derivative wrt e₁e₃ in self.second_der[4]
  • second covariant derivative wrt e₂e₃ in self.second_der[5]
Parameters
  • lmax (int or None, optional): Maximum multipole used in the computation of the derivatives. Default: 3*nside - 1
class SO3Patch(SO3DataField):
554class SO3Patch(SO3DataField):
555    """Class for spin fields in the SO(3) formalism, with Q and U as 2D `np.array` patches. This assumes a flat sky approximation on the patch.
556
557    Parameters
558    ----------
559    Q : np.array
560        Values of the Q field as a 2D array. The x, y dimensions correspond to the ϕ, θ coordinates.
561    
562    U : np.array
563        Values of the U field as a 2D array. Must have the same shape as `Q`.
564        
565    normalise : bool, optional
566        If `True`, the map is normalised to have unit variance.
567        Default: `True`.
568
569    mask : np.array or None, optional
570        Mask where the field is considered. It is a bool array of the same shape that `Q` and `U`.
571        Default: all data is included.
572
573    spacing : float or np.array, optional
574        Spacing between pixels (centres) in each dimension, in radians. If a float, the spacing is the same in all dimensions. 
575        If an array, it must have the same length as the number of dimensions of the field (2).
576        
577    Attributes
578    ----------
579    field : np.array
580        Array of shape `(2, Q.shape)` containing the Q and U values. It can be called as `field(psi)` to obtain the value of $f$ for the given `psi`.
581            
582    dim : int
583        Dimension of the space where the field is defined. In this case, the space is SO(3) and this is 3.
584        
585    name : str
586        Name of the field. In this case, it is `"SO(3) patch"`.
587    
588    first_der : list or None
589        First **covariant** derivatives of the field in an orthonormal basis of the space. It a list of size `3`, each entry has the same structure as `field`.
590    
591    second_der : list or None
592        Second **covariant** derivatives of the field in an orthonormal basis of the space. It a list of size `6`, each entry has the same structure as `field`.
593        The order of the derivatives is diagonal first, i.e., `11`, `22`, `33`, `12`, `13`, `23`.
594        
595    mask : np.array
596        Mask where the field is considered. It is a bool array of the same shape that `Q` and `U`.
597
598    spacing : float or np.array
599        Spacing between pixels (centres) in each dimension. If a float, the spacing is the same in all dimensions.
600        If an array, it must have the same length as the number of dimensions of the field (2).
601
602    """   
603    def __init__(self, Q, U, normalise=True, mask=None, spacing=np.deg2rad(1./60.)):
604        super().__init__(Q, U, name="SO(3) HEALPix map", mask=mask)
605        if self.mask.shape != Q.shape:
606            raise ValueError('The map and the mask have different shapes')
607        
608        if normalise:
609            σ2 = self.get_variance()
610            self.field /= np.sqrt(σ2)
611        self.spacing = spacing
612           
613    def get_variance(self):
614        """Compute the variance of the SO(3) Healpix map within the sky mask. 
615
616        Returns
617        -------
618        var : float
619            The variance of the map within the mask.
620
621        """    
622        return (np.mean(self.field[:,self.mask]**2.))
623    
624    def get_first_der(self):
625        """Compute the covariant first derivatives of the SO(3) Healpix map. Flat sky is assumed.
626        It stores:
627        
628        - first covariant derivative wrt ϕ in self.first_der[0]
629        - first covariant derivative wrt θ in self.first_der[1]
630        - first covariant derivative wrt ψ in self.first_der[2]
631        
632        """    
633        Q_grad = np.gradient(self.field[0], self.spacing, edge_order=2)[::-1]  # order θ,ϕ
634        U_grad = np.gradient(self.field[1], self.spacing, edge_order=2)[::-1]
635
636        fphi = QUarray(Q_grad[1], U_grad[1])
637        ftheta = QUarray(Q_grad[0], U_grad[0])
638        fpsi = self.field.derpsi()
639
640        self.first_der = [ fphi / np.sqrt(2.),                  # ∂e₁  (ϕ)
641                           ftheta / np.sqrt(2.),                # ∂e₂  (θ)
642                           fpsi / np.sqrt(2.)]                  # ∂e₃  (ψ)
643
644        
645            
646    def get_second_der(self):
647        """Compute the covariant second derivatives of the input Healpix scalar map. 
648        It stores:
649        
650        - second covariant derivative wrt ϕϕ in self.second_der[0]
651        - second covariant derivative wrt θθ in self.second_der[1]
652        - second covariant derivative wrt ψψ in self.second_der[2]
653        - second covariant derivative wrt ϕθ in self.second_der[3]
654        - second covariant derivative wrt ϕψ in self.second_der[4]
655        - second covariant derivative wrt θψ in self.second_der[5]
656
657        """    
658        if self.first_der is None:
659            self.get_first_der()
660
661        fphi = self.first_der[0] * np.sqrt(2.)
662        ftheta = self.first_der[1] * np.sqrt(2.)
663        fpsi = self.field.derpsi()
664
665        Q_phi_der = np.gradient(fphi[0], self.spacing, edge_order=2)[::-1]  # order θ,ϕ
666        U_phi_der = np.gradient(fphi[1], self.spacing, edge_order=2)[::-1]  # order θ,ϕ
667        Q_the_der = np.gradient(ftheta[0], self.spacing, edge_order=2)[::-1]
668        U_the_der = np.gradient(ftheta[1], self.spacing, edge_order=2)[::-1]
669
670        fthetatheta = QUarray(Q_the_der[0], U_the_der[0])
671        fthetaphi = QUarray(Q_the_der[1], U_the_der[1])
672        fphiphi = QUarray(Q_phi_der[1], U_phi_der[1])
673        del Q_the_der, U_the_der, Q_phi_der, U_phi_der
674
675        # order 11, 22, 33, 12, 13, 23
676        self.second_der = [ fphiphi / 2.,
677                            fthetatheta / 2.,
678                            fpsi.derpsi() / 2.,
679                            fthetaphi / 2.,
680                            fphi.derpsi() / 2.,
681                            ftheta.derpsi() / 2. ]

Class for spin fields in the SO(3) formalism, with Q and U as 2D np.array patches. This assumes a flat sky approximation on the patch.

Parameters
  • Q (np.array): Values of the Q field as a 2D array. The x, y dimensions correspond to the ϕ, θ coordinates.
  • U (np.array): Values of the U field as a 2D array. Must have the same shape as Q.
  • normalise (bool, optional): If True, the map is normalised to have unit variance. Default: True.
  • mask (np.array or None, optional): Mask where the field is considered. It is a bool array of the same shape that Q and U. Default: all data is included.
  • spacing (float or np.array, optional): Spacing between pixels (centres) in each dimension, in radians. If a float, the spacing is the same in all dimensions. If an array, it must have the same length as the number of dimensions of the field (2).
Attributes
  • field (np.array): Array of shape (2, Q.shape) containing the Q and U values. It can be called as field(psi) to obtain the value of $f$ for the given psi.
  • dim (int): Dimension of the space where the field is defined. In this case, the space is SO(3) and this is 3.
  • name (str): Name of the field. In this case, it is "SO(3) patch".
  • first_der (list or None): First covariant derivatives of the field in an orthonormal basis of the space. It a list of size 3, each entry has the same structure as field.
  • second_der (list or None): Second covariant derivatives of the field in an orthonormal basis of the space. It a list of size 6, each entry has the same structure as field. The order of the derivatives is diagonal first, i.e., 11, 22, 33, 12, 13, 23.
  • mask (np.array): Mask where the field is considered. It is a bool array of the same shape that Q and U.
  • spacing (float or np.array): Spacing between pixels (centres) in each dimension. If a float, the spacing is the same in all dimensions. If an array, it must have the same length as the number of dimensions of the field (2).
SO3Patch(Q, U, normalise=True, mask=None, spacing=0.0002908882086657216)
603    def __init__(self, Q, U, normalise=True, mask=None, spacing=np.deg2rad(1./60.)):
604        super().__init__(Q, U, name="SO(3) HEALPix map", mask=mask)
605        if self.mask.shape != Q.shape:
606            raise ValueError('The map and the mask have different shapes')
607        
608        if normalise:
609            σ2 = self.get_variance()
610            self.field /= np.sqrt(σ2)
611        self.spacing = spacing
def get_variance(self):
613    def get_variance(self):
614        """Compute the variance of the SO(3) Healpix map within the sky mask. 
615
616        Returns
617        -------
618        var : float
619            The variance of the map within the mask.
620
621        """    
622        return (np.mean(self.field[:,self.mask]**2.))

Compute the variance of the SO(3) Healpix map within the sky mask.

Returns
  • var (float): The variance of the map within the mask.
def get_first_der(self):
624    def get_first_der(self):
625        """Compute the covariant first derivatives of the SO(3) Healpix map. Flat sky is assumed.
626        It stores:
627        
628        - first covariant derivative wrt ϕ in self.first_der[0]
629        - first covariant derivative wrt θ in self.first_der[1]
630        - first covariant derivative wrt ψ in self.first_der[2]
631        
632        """    
633        Q_grad = np.gradient(self.field[0], self.spacing, edge_order=2)[::-1]  # order θ,ϕ
634        U_grad = np.gradient(self.field[1], self.spacing, edge_order=2)[::-1]
635
636        fphi = QUarray(Q_grad[1], U_grad[1])
637        ftheta = QUarray(Q_grad[0], U_grad[0])
638        fpsi = self.field.derpsi()
639
640        self.first_der = [ fphi / np.sqrt(2.),                  # ∂e₁  (ϕ)
641                           ftheta / np.sqrt(2.),                # ∂e₂  (θ)
642                           fpsi / np.sqrt(2.)]                  # ∂e₃  (ψ)

Compute the covariant first derivatives of the SO(3) Healpix map. Flat sky is assumed. It stores:

  • first covariant derivative wrt ϕ in self.first_der[0]
  • first covariant derivative wrt θ in self.first_der[1]
  • first covariant derivative wrt ψ in self.first_der[2]
def get_second_der(self):
646    def get_second_der(self):
647        """Compute the covariant second derivatives of the input Healpix scalar map. 
648        It stores:
649        
650        - second covariant derivative wrt ϕϕ in self.second_der[0]
651        - second covariant derivative wrt θθ in self.second_der[1]
652        - second covariant derivative wrt ψψ in self.second_der[2]
653        - second covariant derivative wrt ϕθ in self.second_der[3]
654        - second covariant derivative wrt ϕψ in self.second_der[4]
655        - second covariant derivative wrt θψ in self.second_der[5]
656
657        """    
658        if self.first_der is None:
659            self.get_first_der()
660
661        fphi = self.first_der[0] * np.sqrt(2.)
662        ftheta = self.first_der[1] * np.sqrt(2.)
663        fpsi = self.field.derpsi()
664
665        Q_phi_der = np.gradient(fphi[0], self.spacing, edge_order=2)[::-1]  # order θ,ϕ
666        U_phi_der = np.gradient(fphi[1], self.spacing, edge_order=2)[::-1]  # order θ,ϕ
667        Q_the_der = np.gradient(ftheta[0], self.spacing, edge_order=2)[::-1]
668        U_the_der = np.gradient(ftheta[1], self.spacing, edge_order=2)[::-1]
669
670        fthetatheta = QUarray(Q_the_der[0], U_the_der[0])
671        fthetaphi = QUarray(Q_the_der[1], U_the_der[1])
672        fphiphi = QUarray(Q_phi_der[1], U_phi_der[1])
673        del Q_the_der, U_the_der, Q_phi_der, U_phi_der
674
675        # order 11, 22, 33, 12, 13, 23
676        self.second_der = [ fphiphi / 2.,
677                            fthetatheta / 2.,
678                            fpsi.derpsi() / 2.,
679                            fthetaphi / 2.,
680                            fphi.derpsi() / 2.,
681                            ftheta.derpsi() / 2. ]

Compute the covariant second derivatives of the input Healpix scalar map. It stores:

  • second covariant derivative wrt ϕϕ in self.second_der[0]
  • second covariant derivative wrt θθ in self.second_der[1]
  • second covariant derivative wrt ψψ in self.second_der[2]
  • second covariant derivative wrt ϕθ in self.second_der[3]
  • second covariant derivative wrt ϕψ in self.second_der[4]
  • second covariant derivative wrt θψ in self.second_der[5]
class QUarray(numpy.ndarray):
19class QUarray(np.ndarray):
20    '''Array to store Q and U maps in the SO(3) formalism.
21
22    Parameters
23    ----------
24    Q : np.array
25        Values of the Q field. It can be an array of arrays to represent several maps.
26        
27    U : np.array
28        Values of the U field. Must have the same shape as `Q`.
29        
30    Notes
31    -----
32    The class is a subclass of `np.ndarray`, so it can be used as a numpy array.
33    CAREFUL: the multiplication of two SO3arrays is not the same as the multiplication of the $f$ that they represent. 
34    In order to multiply two SO3arrays, you have to use the __call__ method, i.e. multiply the two SO3arrays after evaluating them at a given positions psi.
35    however, you can multiply a SO3array with a scalar, or add or substract SO3arrays, and the result will be the correct SO3array
36    '''
37    def __new__(cls, Q, U):
38        obj = np.asarray([Q,U]).view(cls)
39        return obj
40
41    def __array_finalize__(self, obj):
42        if obj is None: return
43
44    def __call__(self, psi):
45        '''Return the value of f at a given position psi.
46        
47        Parameters
48        ----------
49        psi : float or np.array
50            The position at which to evaluate f, in radians. It must be broadcastable with the shape of `Q` and `U`.
51            
52        Returns
53        -------
54        f : np.array
55            The value of f at the given position psi.
56        '''
57        return np.asarray(self[0]*np.cos(2.*psi) - self[1]*np.sin(2.*psi))
58    
59    def derpsi(self):
60        '''Return the derivative of f with respect to psi, which can be exactly computed in the SO(3) formalism.
61
62        Returns
63        -------
64        df : np.array
65            The derivative of f with respect to psi.
66        '''
67        return QUarray(-2.*self[1], 2.*self[0])
68    
69    def modulus(self):
70        '''Return the modulus of f.
71
72        Returns
73        -------
74        fmod : np.array
75            The modulus of f.
76        '''
77        return np.sqrt(self[0]**2 + self[1]**2)
78    
79    def pol_angle(self):
80        '''Return the polarization angle of f.
81
82        Returns
83        -------
84        pol_angle : np.array
85            The polarization angle of f.
86        '''
87        angle = np.arctan2(-self[1], self[0])/2
88        angle = angle + np.pi/2*(1-np.sign(self(angle)))/2  # I don't remember why we did this, but it seems to work
89        angle[angle<0] = angle[angle<0] + np.pi
90        return angle

Array to store Q and U maps in the SO(3) formalism.

Parameters
  • Q (np.array): Values of the Q field. It can be an array of arrays to represent several maps.
  • U (np.array): Values of the U field. Must have the same shape as Q.
Notes

The class is a subclass of np.ndarray, so it can be used as a numpy array. CAREFUL: the multiplication of two SO3arrays is not the same as the multiplication of the $f$ that they represent. In order to multiply two SO3arrays, you have to use the __call__ method, i.e. multiply the two SO3arrays after evaluating them at a given positions psi. however, you can multiply a SO3array with a scalar, or add or substract SO3arrays, and the result will be the correct SO3array

def derpsi(self):
59    def derpsi(self):
60        '''Return the derivative of f with respect to psi, which can be exactly computed in the SO(3) formalism.
61
62        Returns
63        -------
64        df : np.array
65            The derivative of f with respect to psi.
66        '''
67        return QUarray(-2.*self[1], 2.*self[0])

Return the derivative of f with respect to psi, which can be exactly computed in the SO(3) formalism.

Returns
  • df (np.array): The derivative of f with respect to psi.
def modulus(self):
69    def modulus(self):
70        '''Return the modulus of f.
71
72        Returns
73        -------
74        fmod : np.array
75            The modulus of f.
76        '''
77        return np.sqrt(self[0]**2 + self[1]**2)

Return the modulus of f.

Returns
  • fmod (np.array): The modulus of f.
def pol_angle(self):
79    def pol_angle(self):
80        '''Return the polarization angle of f.
81
82        Returns
83        -------
84        pol_angle : np.array
85            The polarization angle of f.
86        '''
87        angle = np.arctan2(-self[1], self[0])/2
88        angle = angle + np.pi/2*(1-np.sign(self(angle)))/2  # I don't remember why we did this, but it seems to work
89        angle[angle<0] = angle[angle<0] + np.pi
90        return angle

Return the polarization angle of f.

Returns
  • pol_angle (np.array): The polarization angle of f.
Inherited Members
numpy.ndarray
dumps
dump
all
any
argmax
argmin
argpartition
argsort
astype
byteswap
choose
clip
compress
conj
conjugate
copy
cumprod
cumsum
diagonal
dot
fill
flatten
getfield
item
itemset
max
mean
min
newbyteorder
nonzero
partition
prod
ptp
put
ravel
repeat
reshape
resize
round
searchsorted
setfield
setflags
sort
squeeze
std
sum
swapaxes
take
tobytes
tofile
tolist
tostring
trace
transpose
var
view
ndim
flags
shape
strides
data
itemsize
size
nbytes
base
dtype
real
imag
flat
ctypes
T
class SO3DataField(pynkowski.data.base_da.DataField):
 93class SO3DataField(DataField):
 94    """Class for spherical spin fields in the SO(3) formalism, with Q and U in HEALPix format.
 95
 96    Parameters
 97    ----------
 98    Q : np.array
 99        Values of the Q field.
100    
101    U : np.array
102        Values of the U field. Must have the same shape as `Q`.
103        
104    name : str, optional
105        Name of the field.
106        Defaul : `'SO(3) DataField'`
107    
108    first_der : np.array or None, optional
109        First **covariant** derivatives of the field in an orthonormal basis of the space. Same structure as `field`, and shape `(3, field.shape)`.
110    
111    second_der : np.array or None, optional
112        Second **covariant** derivatives of the field in an orthonormal basis of the space. Same structure as `field`, and shape `(6, field.shape)`.
113        The order of the derivatives is diagonal first, e.g. in `dim=3`: `11`, `22`, `33`, `12`, `13`, `23`.
114        
115    mask : np.array or None, optional
116        Mask where the field is considered. It is a bool array of the same shape that `Q` and `U`.
117        Default: all data is included.
118        
119    Attributes
120    ----------
121    field : np.array
122        Array of shape `(2, Q.shape)` containing the Q and U values. It can be called as `field(psi)` to obtain the value of $f$ for the given `psi`.
123        
124    dim : int
125        Dimension of the space where the field is defined. In this case, the space is SO(3) and this is 3.
126        
127    name : str
128        Name of the field.
129    
130    first_der : list or None
131        First **covariant** derivatives of the field in an orthonormal basis of the space. It a list of size `3`, each entry has the same structure as `field`.
132    
133    second_der : list or None
134        Second **covariant** derivatives of the field in an orthonormal basis of the space. It a list of size `6`, each entry has the same structure as `field`.
135        The order of the derivatives is diagonal first, i.e., `11`, `22`, `33`, `12`, `13`, `23`.
136        
137    mask : np.array
138        Sky mask where the field is considered. It is a bool array of the same shape that `Q` and `U`.
139
140    """   
141    def __init__(self, Q, U, name="SO(3) DataField", first_der=None, second_der=None, mask=None):
142        assert Q.shape == U.shape, "`Q` and `U` must have the same shape"
143        field = QUarray(Q, U)
144        super().__init__(field, 3, name=name, mask=mask, first_der=first_der, second_der=second_der)
145        if mask is None:
146            self.mask = np.ones(Q.shape, dtype=bool)
147        else:
148            assert mask.shape == Q.shape, "`mask` must have the same shape as `Q` and `U`"
149        
150    def get_first_der(self):
151        raise NotImplementedError("First derivatives are not implemented for arbitrary SO(3) DataFields, please use `SO3Healpix` or `SO3Patch`.")
152    
153    def get_second_der(self):
154        raise NotImplementedError("Second derivatives are not implemented for arbitrary SO(3) DataFields, please use `SO3Healpix` or `SO3Patch`.")
155        
156    def _V0(self, us, dus, verbose=True):
157        """Compute the first Minkowski Functional, $V_0$, normalized by the volume of the space. Internal interface for `pynkowski.V0`.
158
159        Parameters
160        ----------
161        us : np.array
162            The thresholds where $V_0$ is computed.
163
164        dus : np.array
165            The bin sizes of `us`. Ignored.
166
167        verbose : bool, optional
168            If True (default), progress bars are shown for the computations on data.
169        
170        Returns
171        -------
172        V0 : np.array()
173            The values of the first Minkowski Functional at the given thresholds.
174        
175        """
176        stat = np.zeros_like(us)
177        modulus = self.field.modulus()
178
179        for ii in tqdm(np.arange(us.shape[0]), disable=not verbose):
180            lenghts = np.zeros_like(modulus)
181            lenghts[us[ii]<=-modulus] = np.pi
182            mask = (us[ii]>-modulus) & (us[ii]<modulus)
183            lenghts[mask] = np.arccos(us[ii]/modulus[mask])
184            stat[ii] = np.mean(lenghts[self.mask])/np.pi
185        return stat
186    
187    def _V1(self, us, dus, verbose=True):
188        """Compute the second Minkowski Functional, $V_1$, normalized by the volume of the space. Internal interface for `pynkowski.V1`.
189
190        Parameters
191        ----------
192        us : np.array
193            The thresholds where $V_1$ is computed.
194            
195        dus : np.array
196            The bin sizes of `us`. Ignored.
197            
198        verbose : bool, optional
199            If True (default), progress bars are shown for the computations on data.
200        
201        Returns
202        -------
203        V1 : np.array()
204            The values of the second Minkowski Functional at the given thresholds.
205        
206        """
207        if self.first_der is None:
208            self.get_first_der()
209        modulus = self.field.modulus()
210        pol_angle = self.field.pol_angle()
211        # theta = get_theta(self.nside)
212        stat = np.zeros_like(us)
213        
214        for ii in tqdm(np.arange(us.shape[0]), disable=not verbose):
215            integrand_1 = np.zeros_like(modulus)
216            integrand_2 = np.zeros_like(modulus)
217            # thetamask = ~(np.isclose(np.cos(theta),0, atol=1.e-2))
218            mask = (us[ii]>-modulus) & (us[ii]<modulus)   # & thetamask
219            if np.sum(mask) == 0:
220                continue
221            first_der_t = QUarray(*zip(*self.first_der))[:,:,mask]
222
223            psi_1 = pol_angle[mask] + np.arccos(us[ii]/modulus[mask])/2.
224            psi_2 = pol_angle[mask] - np.arccos(us[ii]/modulus[mask])/2.
225            
226            integrand_1[mask] = np.sqrt(np.sum(first_der_t(psi_1)**2., axis=0))
227            integrand_2[mask] = np.sqrt(np.sum(first_der_t(psi_2)**2., axis=0))
228            
229            total_integrand = (integrand_1 + integrand_2)   # /np.mean(thetamask)
230            stat[ii] = np.mean(total_integrand[self.mask])/np.pi
231
232        return _MF_prefactor(self.dim, 1) * stat
233    
234    
235    def _V2(self, us, dus, verbose=True):
236        """Compute the second Minkowski Functional, $V_2$, normalized by the volume of the space. Internal interface for `pynkowski.V2`.
237
238        Parameters
239        ----------
240        us : np.array
241            The thresholds where $V_2$ is computed.
242            
243        dus : np.array
244            The bin sizes of `us`. Ignored.
245            
246        verbose : bool, optional
247            If True (default), progress bars are shown for the computations on data.
248        
249        Returns
250        -------
251        V2 : np.array()
252            The values of the second Minkowski Functional at the given thresholds.
253        
254        """
255        if self.first_der is None:
256            self.get_first_der()
257        if self.second_der is None:
258            self.get_second_der()
259            
260        modulus = self.field.modulus()
261        pol_angle = self.field.pol_angle()
262        # theta = get_theta(self.nside)
263        stat = np.zeros_like(us)
264
265        for ii in tqdm(np.arange(us.shape[0]), disable=not verbose):
266            integrand_1 = np.zeros_like(modulus)
267            integrand_2 = np.zeros_like(modulus)
268            
269            # thetamask = ~(np.isclose(np.cos(theta),0, atol=1.e-2))
270            mask = (us[ii]>-modulus) & (us[ii]<modulus)  #& thetamask
271            if np.sum(mask) == 0:
272                continue
273            # pixs = np.arange(12*self.nside**2)[mask]
274            first_der_t = QUarray(*zip(*self.first_der))[:,:,mask]       # This just "transposes" the array without losing QUarray functionality, so that the first index is the Q/U component, the second is the derivative direction, and the third is the pixel
275            second_der_t = QUarray(*zip(*self.second_der))[:,:,mask]
276            
277            psi_1 = pol_angle[mask] + np.arccos(us[ii]/modulus[mask])/2.
278            psi_2 = pol_angle[mask] - np.arccos(us[ii]/modulus[mask])/2.
279
280            integrand_1[mask] = self._H(first_der_t, second_der_t, psi_1)
281            integrand_2[mask] = self._H(first_der_t, second_der_t, psi_2)
282        
283            total_integrand = (integrand_1 + integrand_2)  #/np.mean(thetamask)
284        
285            stat[ii] = np.mean(total_integrand[self.mask])/np.pi
286
287        return _MF_prefactor(self.dim, 2) * stat
288    
289    def _V3(self, us, dus, verbose=True):
290        """Compute the second Minkowski Functional, $V_3$, normalized by the volume of the space. Internal interface for `pynkowski.V3`.
291
292        Parameters
293        ----------
294        us : np.array
295            The thresholds where $V_3$ is computed.
296            
297        dus : np.array
298            The bin sizes of `us`. Ignored.
299            
300        verbose : bool, optional
301            If True (default), progress bars are shown for the computations on data.
302        
303        Returns
304        -------
305        V3 : np.array()
306            The values of the second Minkowski Functional at the given thresholds.
307        
308        """
309        if self.first_der is None:
310            self.get_first_der()
311        if self.second_der is None:
312            self.get_second_der()
313            
314        modulus = self.field.modulus()
315        pol_angle = self.field.pol_angle()
316        # theta = get_theta(self.nside)
317        stat = np.zeros_like(us)
318
319        for ii in tqdm(np.arange(us.shape[0]), disable=not verbose):
320            integrand_1 = np.zeros_like(modulus)
321            integrand_2 = np.zeros_like(modulus)
322            
323            # thetamask = ~(np.isclose(np.cos(theta),0, atol=1.e-2))
324            mask = (us[ii]>-modulus) & (us[ii]<modulus)  #& thetamask
325            if np.sum(mask) == 0:
326                continue
327            # pixs = np.arange(12*self.nside**2)[mask]
328            first_der_t = QUarray(*zip(*self.first_der))[:,:,mask]
329            second_der_t = QUarray(*zip(*self.second_der))[:,:,mask]
330            
331            psi_1 = pol_angle[mask] + np.arccos(us[ii]/modulus[mask])/2.
332            psi_2 = pol_angle[mask] - np.arccos(us[ii]/modulus[mask])/2.
333
334            integrand_1[mask] = self._K(first_der_t, second_der_t, psi_1)
335            integrand_2[mask] = self._K(first_der_t, second_der_t, psi_2)
336        
337            total_integrand = (integrand_1 + integrand_2)  #/np.mean(thetamask)
338        
339            stat[ii] = np.mean(total_integrand[self.mask])/np.pi
340
341        return _MF_prefactor(self.dim, 3) * stat
342
343    @staticmethod
344    def _H(first_der_t, second_der_t, psi):
345        """Compute the mean curvature (times 2) of the field at a given angle.
346        """
347        first_der_t = first_der_t(psi)
348        second_der_t = second_der_t(psi)
349        hess = np.array([[second_der_t[0], second_der_t[3], second_der_t[4]],
350                        [second_der_t[3], second_der_t[1], second_der_t[5]],
351                        [second_der_t[4], second_der_t[5], second_der_t[2]]])
352        norm_grad = first_der_t / np.sqrt(np.sum(first_der_t**2., axis=0))
353        return np.einsum('j...,jk...,k... -> ...', norm_grad, hess, norm_grad) - np.trace(hess, axis1=0, axis2=1)
354    
355    @staticmethod
356    def _K(first_der_t, second_der_t, psi):
357        """Compute the Gaussian curvature of the field at a given angle.
358        """
359        first_der_t = first_der_t(psi)
360        mod_grad = np.sqrt(np.sum(first_der_t**2., axis=0))
361        second_der_t = second_der_t(psi) / mod_grad
362        norm_grad = first_der_t / mod_grad
363        extended_hessian_det = np.linalg.det(np.array([[second_der_t[0], second_der_t[3], second_der_t[4], norm_grad[0]],
364                                                        [second_der_t[3], second_der_t[1], second_der_t[5], norm_grad[1]],
365                                                        [second_der_t[4], second_der_t[5], second_der_t[2], norm_grad[2]],
366                                                        [norm_grad[0], norm_grad[1], norm_grad[2], np.zeros_like(norm_grad[0])]]).T).T
367        return -extended_hessian_det*mod_grad

Class for spherical spin fields in the SO(3) formalism, with Q and U in HEALPix format.

Parameters
  • Q (np.array): Values of the Q field.
  • U (np.array): Values of the U field. Must have the same shape as Q.
  • name (str, optional): Name of the field. Defaul : 'SO(3) DataField'
  • first_der (np.array or None, optional): First covariant derivatives of the field in an orthonormal basis of the space. Same structure as field, and shape (3, field.shape).
  • second_der (np.array or None, optional): Second covariant derivatives of the field in an orthonormal basis of the space. Same structure as field, and shape (6, field.shape). The order of the derivatives is diagonal first, e.g. in dim=3: 11, 22, 33, 12, 13, 23.
  • mask (np.array or None, optional): Mask where the field is considered. It is a bool array of the same shape that Q and U. Default: all data is included.
Attributes
  • field (np.array): Array of shape (2, Q.shape) containing the Q and U values. It can be called as field(psi) to obtain the value of $f$ for the given psi.
  • dim (int): Dimension of the space where the field is defined. In this case, the space is SO(3) and this is 3.
  • name (str): Name of the field.
  • first_der (list or None): First covariant derivatives of the field in an orthonormal basis of the space. It a list of size 3, each entry has the same structure as field.
  • second_der (list or None): Second covariant derivatives of the field in an orthonormal basis of the space. It a list of size 6, each entry has the same structure as field. The order of the derivatives is diagonal first, i.e., 11, 22, 33, 12, 13, 23.
  • mask (np.array): Sky mask where the field is considered. It is a bool array of the same shape that Q and U.
SO3DataField( Q, U, name='SO(3) DataField', first_der=None, second_der=None, mask=None)
141    def __init__(self, Q, U, name="SO(3) DataField", first_der=None, second_der=None, mask=None):
142        assert Q.shape == U.shape, "`Q` and `U` must have the same shape"
143        field = QUarray(Q, U)
144        super().__init__(field, 3, name=name, mask=mask, first_der=first_der, second_der=second_der)
145        if mask is None:
146            self.mask = np.ones(Q.shape, dtype=bool)
147        else:
148            assert mask.shape == Q.shape, "`mask` must have the same shape as `Q` and `U`"
def get_first_der(self):
150    def get_first_der(self):
151        raise NotImplementedError("First derivatives are not implemented for arbitrary SO(3) DataFields, please use `SO3Healpix` or `SO3Patch`.")
def get_second_der(self):
153    def get_second_der(self):
154        raise NotImplementedError("Second derivatives are not implemented for arbitrary SO(3) DataFields, please use `SO3Healpix` or `SO3Patch`.")