pynkowski.theory.spingaussian
This submodule contains the definition for spin field with isotropic Gaussian Q and U.
1"""This submodule contains the definition for spin field with isotropic Gaussian Q and U. 2""" 3import numpy as np 4from .utils_th import get_μ, get_σ, lkc_ambient_dict, get_C2, flag, rho 5from .gaussian import Gaussian 6 7from .base_th import _prepare_lkc 8 9 10 11def LKC_spin(j, us, mu, dim=3, lkc_ambient=lkc_ambient_dict["SO3"], Ks=[1., 0.31912, 0.7088, 1.]): 12 """Compute the expected value of the Lipschitz–Killing Curvatures (LKC) of the excursion set for Gaussian Isotropic fields. 13 14 Parameters 15 ---------- 16 j : int 17 The index of the LKC. 18 19 us : np.array 20 The thresholds where the LKC is evaluated. 21 22 mu : float 23 The value of the derivative of the covariance function at the origin for the U and Q fields. 24 25 dim : int, optional 26 The dimension of the ambient manifold, SO(3). Must be 3. 27 28 lkc_ambient : np.array, optional 29 An array of the Lipschitz–Killing Curvatures of the ambient manifold. Its lenght must be `dim+1` if `dim` is also given. 30 31 Ks : np.array, optional 32 The normalisation factors for the Minkowski Functionals. 33 34 Returns 35 ---------- 36 LKC : np.array 37 The expected value of the Lipschitz–Killing Curvatures at the thresholds. 38 """ 39 dim, lkc_ambient = _prepare_lkc(dim, lkc_ambient) 40 assert dim==3, "The LKC for spin fields is only defined on SO(3), which has `dim=3`" 41 result = np.zeros_like(us) 42 KLs = Ks[::-1] # To get the order as the LKC index, not the MF index 43 KLs[-1] *= 2.**1.5 / 5. # This is the determinant of the covariance matrix of the derivatives of the field 44 KLs /= np.array([1., 1., mu**0.5, mu]) 45 for k in np.arange(0,dim-j+1): 46 result += flag(k+j, k) * rho(k, us) * lkc_ambient[k+j] / KLs[k+j] 47 return result * KLs[j] /lkc_ambient[-1] 48 49 50 51 52class SpinGaussian(Gaussian): 53 """Class for Spin Isotropic Gaussian fields in the SO(3) formalism. 54 55 Parameters 56 ---------- 57 cls : np.array 58 Angular power spectrum of the field (cls_E + cls_B). 59 60 normalise : bool, optional 61 If `True`, normalise the field to unit variance. 62 Default : True 63 64 fsky : float, optional 65 Fraction of the sky covered by the field, `0<fsky<=1`. 66 Default : 1. 67 68 Ks : np.array, optional 69 The normalisation constants for the MFs of the field: [K_0, K_1, K_2, K_3]. 70 Default : [1., 0.31912, 0.7088, 1.] as found in Carrón Duque et al (2023) 71 72 leading_order : bool, optional 73 Whether to use only the leading order in μ for the computation of the MFs or the exact expression (with two terms). 74 Default : False (exact expression) 75 76 Attributes 77 ---------- 78 cls : np.array 79 Angular Power Spectrum of the field (cls_E + cls_B). 80 81 fsky : float 82 Fraction of the sky covered by the field. 83 84 dim : int 85 Dimension of the space where the field is defined, in this case this is 3. 86 87 name : str 88 Name of the field, `"Spin Isotropic Gaussian"` by default. 89 90 sigma : float 91 The standard deviation of the field. 92 93 mu : float 94 The derivative of the covariance function at the origin, times $-2$ (in the spatial coordinates θϕ). Equal to the variance of the first derivatives of the field. 95 96 nu : float 97 The second derivative of the covariance function at the origin (in the spatial coordinates θϕ). 98 99 C2 : float 100 The second derivative of the angular covariance function at 1 (in the spatial coordinates θϕ). 101 102 lkc_ambient : np.array or None 103 The values for the Lipschitz–Killing Curvatures of the ambient space. 104 105 Ks : np.array 106 The normalisation constants for the MFs of the field: [K_0, K_1, K_2, K_3]. 107 108 leading_order : bool 109 Whether to use only the leading order in μ for the computation of the MFs or the exact expression (with two terms). 110 111 """ 112 def __init__(self, cls, normalise=True, fsky=1., Ks=[1., 0.31912, 0.7088, 1.], leading_order=True): 113 if normalise: 114 cls /= get_σ(cls) 115 self.sigma = 1. 116 else: 117 self.sigma = np.sqrt(get_σ(cls)) 118 self.cls = cls 119 self.fsky = fsky 120 self.mu = get_μ(cls) 121 self.C2 = get_C2(cls) 122 self.nu = self.C2/4. - self.mu/24. 123 super().__init__(3, sigma=self.sigma, mu=self.mu, nu=self.nu, lkc_ambient=lkc_ambient_dict["SO3"]*self.fsky) 124 self.leading_order = leading_order 125 if leading_order: 126 self.lkc_ambient[1] = 0. 127 self.name = 'Spin Isotropic Gaussian' 128 self.Ks = Ks 129 130 def LKC(self, j, us): 131 """Compute the expected values of the Lipschitz–Killing Curvatures of the excursion sets at thresholds `us`, $\mathbb{L}_j(A_u(f))$. 132 133 Parameters 134 ---------- 135 j : int 136 Index of the LKC to compute, `0 < j < dim`. 137 138 us : np.array 139 The thresholds considered for the computation. 140 141 Returns 142 ------- 143 lkc : np.array 144 Expected value of LKC evaluated at the thresholds. 145 146 """ 147 us /= self.sigma 148 return LKC_spin(j, us, self.mu, dim=self.dim, lkc_ambient=self.lkc_ambient, Ks=self.Ks) 149 150 151__all__ = ["SpinGaussian"]
53class SpinGaussian(Gaussian): 54 """Class for Spin Isotropic Gaussian fields in the SO(3) formalism. 55 56 Parameters 57 ---------- 58 cls : np.array 59 Angular power spectrum of the field (cls_E + cls_B). 60 61 normalise : bool, optional 62 If `True`, normalise the field to unit variance. 63 Default : True 64 65 fsky : float, optional 66 Fraction of the sky covered by the field, `0<fsky<=1`. 67 Default : 1. 68 69 Ks : np.array, optional 70 The normalisation constants for the MFs of the field: [K_0, K_1, K_2, K_3]. 71 Default : [1., 0.31912, 0.7088, 1.] as found in Carrón Duque et al (2023) 72 73 leading_order : bool, optional 74 Whether to use only the leading order in μ for the computation of the MFs or the exact expression (with two terms). 75 Default : False (exact expression) 76 77 Attributes 78 ---------- 79 cls : np.array 80 Angular Power Spectrum of the field (cls_E + cls_B). 81 82 fsky : float 83 Fraction of the sky covered by the field. 84 85 dim : int 86 Dimension of the space where the field is defined, in this case this is 3. 87 88 name : str 89 Name of the field, `"Spin Isotropic Gaussian"` by default. 90 91 sigma : float 92 The standard deviation of the field. 93 94 mu : float 95 The derivative of the covariance function at the origin, times $-2$ (in the spatial coordinates θϕ). Equal to the variance of the first derivatives of the field. 96 97 nu : float 98 The second derivative of the covariance function at the origin (in the spatial coordinates θϕ). 99 100 C2 : float 101 The second derivative of the angular covariance function at 1 (in the spatial coordinates θϕ). 102 103 lkc_ambient : np.array or None 104 The values for the Lipschitz–Killing Curvatures of the ambient space. 105 106 Ks : np.array 107 The normalisation constants for the MFs of the field: [K_0, K_1, K_2, K_3]. 108 109 leading_order : bool 110 Whether to use only the leading order in μ for the computation of the MFs or the exact expression (with two terms). 111 112 """ 113 def __init__(self, cls, normalise=True, fsky=1., Ks=[1., 0.31912, 0.7088, 1.], leading_order=True): 114 if normalise: 115 cls /= get_σ(cls) 116 self.sigma = 1. 117 else: 118 self.sigma = np.sqrt(get_σ(cls)) 119 self.cls = cls 120 self.fsky = fsky 121 self.mu = get_μ(cls) 122 self.C2 = get_C2(cls) 123 self.nu = self.C2/4. - self.mu/24. 124 super().__init__(3, sigma=self.sigma, mu=self.mu, nu=self.nu, lkc_ambient=lkc_ambient_dict["SO3"]*self.fsky) 125 self.leading_order = leading_order 126 if leading_order: 127 self.lkc_ambient[1] = 0. 128 self.name = 'Spin Isotropic Gaussian' 129 self.Ks = Ks 130 131 def LKC(self, j, us): 132 """Compute the expected values of the Lipschitz–Killing Curvatures of the excursion sets at thresholds `us`, $\mathbb{L}_j(A_u(f))$. 133 134 Parameters 135 ---------- 136 j : int 137 Index of the LKC to compute, `0 < j < dim`. 138 139 us : np.array 140 The thresholds considered for the computation. 141 142 Returns 143 ------- 144 lkc : np.array 145 Expected value of LKC evaluated at the thresholds. 146 147 """ 148 us /= self.sigma 149 return LKC_spin(j, us, self.mu, dim=self.dim, lkc_ambient=self.lkc_ambient, Ks=self.Ks)
Class for Spin Isotropic Gaussian fields in the SO(3) formalism.
Parameters
cls : np.array Angular power spectrum of the field (cls_E + cls_B).
normalise : bool, optional
If True
, normalise the field to unit variance.
Default : True
fsky : float, optional
Fraction of the sky covered by the field, 0<fsky<=1
.
Default : 1.
Ks : np.array, optional The normalisation constants for the MFs of the field: [K_0, K_1, K_2, K_3]. Default : [1., 0.31912, 0.7088, 1.] as found in Carrón Duque et al (2023)
leading_order : bool, optional Whether to use only the leading order in μ for the computation of the MFs or the exact expression (with two terms). Default : False (exact expression)
Attributes
cls : np.array Angular Power Spectrum of the field (cls_E + cls_B).
fsky : float Fraction of the sky covered by the field.
dim : int Dimension of the space where the field is defined, in this case this is 3.
name : str
Name of the field, "Spin Isotropic Gaussian"
by default.
sigma : float The standard deviation of the field.
mu : float The derivative of the covariance function at the origin, times $-2$ (in the spatial coordinates θϕ). Equal to the variance of the first derivatives of the field.
nu : float The second derivative of the covariance function at the origin (in the spatial coordinates θϕ).
C2 : float The second derivative of the angular covariance function at 1 (in the spatial coordinates θϕ).
lkc_ambient : np.array or None The values for the Lipschitz–Killing Curvatures of the ambient space.
Ks : np.array The normalisation constants for the MFs of the field: [K_0, K_1, K_2, K_3].
leading_order : bool Whether to use only the leading order in μ for the computation of the MFs or the exact expression (with two terms).
113 def __init__(self, cls, normalise=True, fsky=1., Ks=[1., 0.31912, 0.7088, 1.], leading_order=True): 114 if normalise: 115 cls /= get_σ(cls) 116 self.sigma = 1. 117 else: 118 self.sigma = np.sqrt(get_σ(cls)) 119 self.cls = cls 120 self.fsky = fsky 121 self.mu = get_μ(cls) 122 self.C2 = get_C2(cls) 123 self.nu = self.C2/4. - self.mu/24. 124 super().__init__(3, sigma=self.sigma, mu=self.mu, nu=self.nu, lkc_ambient=lkc_ambient_dict["SO3"]*self.fsky) 125 self.leading_order = leading_order 126 if leading_order: 127 self.lkc_ambient[1] = 0. 128 self.name = 'Spin Isotropic Gaussian' 129 self.Ks = Ks
131 def LKC(self, j, us): 132 """Compute the expected values of the Lipschitz–Killing Curvatures of the excursion sets at thresholds `us`, $\mathbb{L}_j(A_u(f))$. 133 134 Parameters 135 ---------- 136 j : int 137 Index of the LKC to compute, `0 < j < dim`. 138 139 us : np.array 140 The thresholds considered for the computation. 141 142 Returns 143 ------- 144 lkc : np.array 145 Expected value of LKC evaluated at the thresholds. 146 147 """ 148 us /= self.sigma 149 return LKC_spin(j, us, self.mu, dim=self.dim, lkc_ambient=self.lkc_ambient, Ks=self.Ks)
Compute the expected values of the Lipschitz–Killing Curvatures of the excursion sets at thresholds us
, $\mathbb{L}_j(A_u(f))$.
Parameters
j : int
Index of the LKC to compute, 0 < j < dim
.
us : np.array The thresholds considered for the computation.
Returns
lkc : np.array Expected value of LKC evaluated at the thresholds.