pynkowski.theory.chi2

This submodule contains the classes for the theoretical $\chi^2$ fields.

  1'''This submodule contains the classes for the theoretical $\chi^2$ fields.'''
  2import numpy as np
  3from scipy.special import comb
  4from .utils_th import LKC_Chi2, get_σ, get_μ, lkc_ambient_dict
  5from .base_th import TheoryField
  6        
  7class Chi2(TheoryField):
  8    """General class for Isotropic Chi² fields, to be used directly or as the base for specific Chi² fields.
  9
 10    Parameters
 11    ----------
 12    dim : int
 13        Dimension of the space where the field is defined.
 14        
 15    dof : int, optional
 16        Number of degrees of freedom $dof$ of the field: $\Chi^2 = \sum_{i=0}{dof} N_i^2$.
 17        Default: 2
 18    
 19    sigma : float, optional
 20        The standard deviation of the constituent Gaussian fields $N_i$, such that $\Chi^2 = \sum_{i=0}{dof} N_i^2$.
 21        Default: 1.
 22        
 23    mu : float, optional
 24        The derivative of the covariance function at the origin for the constituent Gaussian fields.
 25        Default: 1.
 26
 27    nu : float, optional
 28        The second derivative of the covariance function at the origin for the constituent Gaussian fields.
 29        Default: 1.
 30        
 31    lkc_ambient : np.array or None, optional
 32        The values for the Lipschitz–Killing Curvatures of the ambient space. If `None`, it is assumed that the volume is 1 and the rest is 0. This is exact for many spaces like Euclidean spaces or the sphere, and exact to leading order in μ for the rest.
 33        Default: None
 34        
 35    Attributes
 36    ----------
 37    dim : int
 38        Dimension of the space where the field is defined.
 39    
 40    name : str
 41        Name of the field, `"Isotropic Chi²"` by default.
 42
 43    dof : int, optional
 44        Number of degrees of freedom $dof$ of the field: $\Chi^2 = \sum_{i=0}{dof} N_i^2$.
 45    
 46    sigma : float
 47        The standard deviation of the field.
 48        
 49    mu : float
 50        The derivative of the covariance function at the origin.
 51
 52    nu : float, optional
 53        The second derivative of the covariance function at the origin for the constituent Gaussian fields.
 54        
 55    lkc_ambient : np.array or None
 56        The values for the Lipschitz–Killing Curvatures of the ambient space.
 57        
 58    """   
 59    def __init__(self, dim, dof=2, sigma=1., mu=1., nu=1., lkc_ambient=None):
 60        super().__init__(dim, name='Isotropic Chi²', sigma=sigma, mu=mu, nu=nu, lkc_ambient=lkc_ambient)
 61        self.dof = dof
 62        
 63    def LKC(self, j, us):
 64        """Compute the expected values of the Lipschitz–Killing Curvatures of the excursion sets at thresholds `us`, $\mathbb{L}_j(A_u(f))$.
 65
 66        Parameters
 67        ----------
 68        j : int
 69            Index of the LKC to compute, `0 < j < dim`.
 70            
 71        us : np.array
 72            The thresholds considered for the computation.
 73
 74        Returns
 75        -------
 76        lkc : np.array
 77            Expected value of LKC evaluated at the thresholds.
 78
 79        """
 80        us /= self.sigma
 81        return LKC_Chi2(j, us, self.mu, dof=self.dof, dim=self.dim, lkc_ambient=self.lkc_ambient)/self.lkc_ambient[-1]
 82    
 83    def MF(self, j, us):
 84        """Compute the expected values of the Minkowski Functionals of the excursion sets at thresholds `us`, $V_j(A_u(f))$.
 85
 86        Parameters
 87        ----------
 88        j : int
 89            Index of the MF to compute, `0 < j < dim`.
 90            
 91        us : np.array
 92            The thresholds considered for the computation.
 93
 94        Returns
 95        -------
 96        mf : np.array
 97            Expected value of MFs evaluated at the thresholds.
 98
 99        """
100        us /= self.sigma
101        return 1./comb(self.dim,j)*LKC_Chi2(self.dim-j, us, self.mu, dof=self.dof, dim=self.dim, lkc_ambient=self.lkc_ambient)/self.lkc_ambient[-1]
102    
103    def V0(self, us):
104        """Compute the expected values of the first Minkowski Functionals of the excursion sets at thresholds `us`, $V_0(A_u(f))$.
105
106        Parameters
107        ----------
108        us : np.array
109            The thresholds considered for the computation.
110
111        Returns
112        -------
113        v0 : np.array
114            Expected value of $V_0$ evaluated at the thresholds.
115
116        """
117        us /= self.sigma
118        return self.MF(0, us)
119    
120    def V1(self, us):
121        """Compute the expected values of the first Minkowski Functionals of the excursion sets at thresholds `us`, $V_1(A_u(f))$.
122
123        Parameters
124        ----------
125        us : np.array
126            The thresholds considered for the computation.
127
128        Returns
129        -------
130        v1 : np.array
131            Expected value of $V_1$ evaluated at the thresholds.
132
133        """
134        us /= self.sigma
135        return self.MF(1, us)
136    
137    def V2(self, us):
138        """Compute the expected values of the first Minkowski Functionals of the excursion sets at thresholds `us`, $V_2(A_u(f))$. Only for fields with `dim ≥ 2`.
139
140        Parameters
141        ----------
142        us : np.array
143            The thresholds considered for the computation.
144
145        Returns
146        -------
147        v2 : np.array
148            Expected value of $V_2$ evaluated at the thresholds.
149
150        """
151        assert self.dim>=2, 'V2 is defined only for fields with dim≥2'
152        us /= self.sigma
153        return self.MF(2, us)
154    
155    def V3(self, us):
156        """Compute the expected values of the first Minkowski Functionals of the excursion sets at thresholds `us`, $V_3(A_u(f))$. Only for fields with `dim ≥ 3`.
157
158        Parameters
159        ----------
160        us : np.array
161            The thresholds considered for the computation.
162
163        Returns
164        -------
165        v3 : np.array
166            Expected value of $V_3$ evaluated at the thresholds.
167
168        """
169        assert self.dim>=3, 'V3 is defined only for fields with dim≥3'
170        us /= self.sigma
171        return self.MF(3, us)
172    
173    # def V4(self, us):
174    #     assert self.dim>=4, 'V4 is defined only for fields with dim≥4'
175    #     us /= self.sigma
176    #     return self.MF(4, us)
177    
178    
179class SphericalChi2(Chi2):
180    """Class for Spherical Isotropic Gaussian fields, to be used directly or as the base for specific Gaussian fields.
181
182    Parameters
183    ----------
184    cls : np.array
185        Angular power spectrum of the polarization field (EE, BB). The shape is `(2,lmax+1)`
186    
187    normalise : bool, optional
188        If `True`, normalise the Q and U fields to unit variance
189        Default : True
190    
191    Attributes
192    ----------
193    cls : np.array
194        Angular Power Spectrum of the polarization field.
195    
196    dim : int
197        Dimension of the space where the field is defined, in this case this is 2.
198    
199    name : str
200        Name of the field, `"Spherical Isotropic Chi²"` by default.
201
202    dof : int, optional
203        Number of degrees of freedom $dof$ of the field, $2$.
204    
205    sigma : float
206        The standard deviation of the field.
207        
208    mu : float
209        The derivative of the covariance function at the origin.
210
211    nu : float, optional
212        The second derivative of the covariance function at the origin for the constituent Gaussian fields.
213        
214    lkc_ambient : np.array or None
215        The values for the Lipschitz–Killing Curvatures of the ambient space.
216        
217    """   
218    def __init__(self, cls, normalise=True):
219        if normalise:
220            cls /= get_σ((cls[0]+cls[1])/2.)
221            self.sigma = 1.
222        else:
223            self.sigma = np.sqrt(get_σ(cls))
224        self.cls = cls
225        self.mu = get_μ(cls[0]+cls[1])
226        super().__init__(dim=2, dof=2, sigma=self.sigma, mu=self.mu, lkc_ambient=lkc_ambient_dict["sphere"])
227        self.name = 'Spherical Isotropic Chi²'        
228        
229
230
231__all__ = ["SphericalChi2", "Chi2"]
232
233__docformat__ = "numpy"
234
235 
class SphericalChi2(Chi2):
180class SphericalChi2(Chi2):
181    """Class for Spherical Isotropic Gaussian fields, to be used directly or as the base for specific Gaussian fields.
182
183    Parameters
184    ----------
185    cls : np.array
186        Angular power spectrum of the polarization field (EE, BB). The shape is `(2,lmax+1)`
187    
188    normalise : bool, optional
189        If `True`, normalise the Q and U fields to unit variance
190        Default : True
191    
192    Attributes
193    ----------
194    cls : np.array
195        Angular Power Spectrum of the polarization field.
196    
197    dim : int
198        Dimension of the space where the field is defined, in this case this is 2.
199    
200    name : str
201        Name of the field, `"Spherical Isotropic Chi²"` by default.
202
203    dof : int, optional
204        Number of degrees of freedom $dof$ of the field, $2$.
205    
206    sigma : float
207        The standard deviation of the field.
208        
209    mu : float
210        The derivative of the covariance function at the origin.
211
212    nu : float, optional
213        The second derivative of the covariance function at the origin for the constituent Gaussian fields.
214        
215    lkc_ambient : np.array or None
216        The values for the Lipschitz–Killing Curvatures of the ambient space.
217        
218    """   
219    def __init__(self, cls, normalise=True):
220        if normalise:
221            cls /= get_σ((cls[0]+cls[1])/2.)
222            self.sigma = 1.
223        else:
224            self.sigma = np.sqrt(get_σ(cls))
225        self.cls = cls
226        self.mu = get_μ(cls[0]+cls[1])
227        super().__init__(dim=2, dof=2, sigma=self.sigma, mu=self.mu, lkc_ambient=lkc_ambient_dict["sphere"])
228        self.name = 'Spherical Isotropic Chi²'        

Class for Spherical Isotropic Gaussian fields, to be used directly or as the base for specific Gaussian fields.

Parameters
  • cls (np.array): Angular power spectrum of the polarization field (EE, BB). The shape is (2,lmax+1)
  • normalise (bool, optional): If True, normalise the Q and U fields to unit variance Default : True
Attributes
  • cls (np.array): Angular Power Spectrum of the polarization field.
  • dim (int): Dimension of the space where the field is defined, in this case this is 2.
  • name (str): Name of the field, "Spherical Isotropic Chi²" by default.
  • dof (int, optional): Number of degrees of freedom $dof$ of the field, $2$.
  • sigma (float): The standard deviation of the field.
  • mu (float): The derivative of the covariance function at the origin.
  • nu (float, optional): The second derivative of the covariance function at the origin for the constituent Gaussian fields.
  • lkc_ambient (np.array or None): The values for the Lipschitz–Killing Curvatures of the ambient space.
SphericalChi2(cls, normalise=True)
219    def __init__(self, cls, normalise=True):
220        if normalise:
221            cls /= get_σ((cls[0]+cls[1])/2.)
222            self.sigma = 1.
223        else:
224            self.sigma = np.sqrt(get_σ(cls))
225        self.cls = cls
226        self.mu = get_μ(cls[0]+cls[1])
227        super().__init__(dim=2, dof=2, sigma=self.sigma, mu=self.mu, lkc_ambient=lkc_ambient_dict["sphere"])
228        self.name = 'Spherical Isotropic Chi²'        
Inherited Members
Chi2
LKC
MF
V0
V1
V2
V3
class Chi2(pynkowski.theory.base_th.TheoryField):
  8class Chi2(TheoryField):
  9    """General class for Isotropic Chi² fields, to be used directly or as the base for specific Chi² fields.
 10
 11    Parameters
 12    ----------
 13    dim : int
 14        Dimension of the space where the field is defined.
 15        
 16    dof : int, optional
 17        Number of degrees of freedom $dof$ of the field: $\Chi^2 = \sum_{i=0}{dof} N_i^2$.
 18        Default: 2
 19    
 20    sigma : float, optional
 21        The standard deviation of the constituent Gaussian fields $N_i$, such that $\Chi^2 = \sum_{i=0}{dof} N_i^2$.
 22        Default: 1.
 23        
 24    mu : float, optional
 25        The derivative of the covariance function at the origin for the constituent Gaussian fields.
 26        Default: 1.
 27
 28    nu : float, optional
 29        The second derivative of the covariance function at the origin for the constituent Gaussian fields.
 30        Default: 1.
 31        
 32    lkc_ambient : np.array or None, optional
 33        The values for the Lipschitz–Killing Curvatures of the ambient space. If `None`, it is assumed that the volume is 1 and the rest is 0. This is exact for many spaces like Euclidean spaces or the sphere, and exact to leading order in μ for the rest.
 34        Default: None
 35        
 36    Attributes
 37    ----------
 38    dim : int
 39        Dimension of the space where the field is defined.
 40    
 41    name : str
 42        Name of the field, `"Isotropic Chi²"` by default.
 43
 44    dof : int, optional
 45        Number of degrees of freedom $dof$ of the field: $\Chi^2 = \sum_{i=0}{dof} N_i^2$.
 46    
 47    sigma : float
 48        The standard deviation of the field.
 49        
 50    mu : float
 51        The derivative of the covariance function at the origin.
 52
 53    nu : float, optional
 54        The second derivative of the covariance function at the origin for the constituent Gaussian fields.
 55        
 56    lkc_ambient : np.array or None
 57        The values for the Lipschitz–Killing Curvatures of the ambient space.
 58        
 59    """   
 60    def __init__(self, dim, dof=2, sigma=1., mu=1., nu=1., lkc_ambient=None):
 61        super().__init__(dim, name='Isotropic Chi²', sigma=sigma, mu=mu, nu=nu, lkc_ambient=lkc_ambient)
 62        self.dof = dof
 63        
 64    def LKC(self, j, us):
 65        """Compute the expected values of the Lipschitz–Killing Curvatures of the excursion sets at thresholds `us`, $\mathbb{L}_j(A_u(f))$.
 66
 67        Parameters
 68        ----------
 69        j : int
 70            Index of the LKC to compute, `0 < j < dim`.
 71            
 72        us : np.array
 73            The thresholds considered for the computation.
 74
 75        Returns
 76        -------
 77        lkc : np.array
 78            Expected value of LKC evaluated at the thresholds.
 79
 80        """
 81        us /= self.sigma
 82        return LKC_Chi2(j, us, self.mu, dof=self.dof, dim=self.dim, lkc_ambient=self.lkc_ambient)/self.lkc_ambient[-1]
 83    
 84    def MF(self, j, us):
 85        """Compute the expected values of the Minkowski Functionals of the excursion sets at thresholds `us`, $V_j(A_u(f))$.
 86
 87        Parameters
 88        ----------
 89        j : int
 90            Index of the MF to compute, `0 < j < dim`.
 91            
 92        us : np.array
 93            The thresholds considered for the computation.
 94
 95        Returns
 96        -------
 97        mf : np.array
 98            Expected value of MFs evaluated at the thresholds.
 99
100        """
101        us /= self.sigma
102        return 1./comb(self.dim,j)*LKC_Chi2(self.dim-j, us, self.mu, dof=self.dof, dim=self.dim, lkc_ambient=self.lkc_ambient)/self.lkc_ambient[-1]
103    
104    def V0(self, us):
105        """Compute the expected values of the first Minkowski Functionals of the excursion sets at thresholds `us`, $V_0(A_u(f))$.
106
107        Parameters
108        ----------
109        us : np.array
110            The thresholds considered for the computation.
111
112        Returns
113        -------
114        v0 : np.array
115            Expected value of $V_0$ evaluated at the thresholds.
116
117        """
118        us /= self.sigma
119        return self.MF(0, us)
120    
121    def V1(self, us):
122        """Compute the expected values of the first Minkowski Functionals of the excursion sets at thresholds `us`, $V_1(A_u(f))$.
123
124        Parameters
125        ----------
126        us : np.array
127            The thresholds considered for the computation.
128
129        Returns
130        -------
131        v1 : np.array
132            Expected value of $V_1$ evaluated at the thresholds.
133
134        """
135        us /= self.sigma
136        return self.MF(1, us)
137    
138    def V2(self, us):
139        """Compute the expected values of the first Minkowski Functionals of the excursion sets at thresholds `us`, $V_2(A_u(f))$. Only for fields with `dim ≥ 2`.
140
141        Parameters
142        ----------
143        us : np.array
144            The thresholds considered for the computation.
145
146        Returns
147        -------
148        v2 : np.array
149            Expected value of $V_2$ evaluated at the thresholds.
150
151        """
152        assert self.dim>=2, 'V2 is defined only for fields with dim≥2'
153        us /= self.sigma
154        return self.MF(2, us)
155    
156    def V3(self, us):
157        """Compute the expected values of the first Minkowski Functionals of the excursion sets at thresholds `us`, $V_3(A_u(f))$. Only for fields with `dim ≥ 3`.
158
159        Parameters
160        ----------
161        us : np.array
162            The thresholds considered for the computation.
163
164        Returns
165        -------
166        v3 : np.array
167            Expected value of $V_3$ evaluated at the thresholds.
168
169        """
170        assert self.dim>=3, 'V3 is defined only for fields with dim≥3'
171        us /= self.sigma
172        return self.MF(3, us)
173    
174    # def V4(self, us):
175    #     assert self.dim>=4, 'V4 is defined only for fields with dim≥4'
176    #     us /= self.sigma
177    #     return self.MF(4, us)

General class for Isotropic Chi² fields, to be used directly or as the base for specific Chi² fields.

Parameters
  • dim (int): Dimension of the space where the field is defined.
  • dof (int, optional): Number of degrees of freedom $dof$ of the field: $\Chi^2 = \sum_{i=0}{dof} N_i^2$. Default: 2
  • sigma (float, optional): The standard deviation of the constituent Gaussian fields $N_i$, such that $\Chi^2 = \sum_{i=0}{dof} N_i^2$. Default: 1.
  • mu (float, optional): The derivative of the covariance function at the origin for the constituent Gaussian fields. Default: 1.
  • nu (float, optional): The second derivative of the covariance function at the origin for the constituent Gaussian fields. Default: 1.
  • lkc_ambient (np.array or None, optional): The values for the Lipschitz–Killing Curvatures of the ambient space. If None, it is assumed that the volume is 1 and the rest is 0. This is exact for many spaces like Euclidean spaces or the sphere, and exact to leading order in μ for the rest. Default: None
Attributes
  • dim (int): Dimension of the space where the field is defined.
  • name (str): Name of the field, "Isotropic Chi²" by default.
  • dof (int, optional): Number of degrees of freedom $dof$ of the field: $\Chi^2 = \sum_{i=0}{dof} N_i^2$.
  • sigma (float): The standard deviation of the field.
  • mu (float): The derivative of the covariance function at the origin.
  • nu (float, optional): The second derivative of the covariance function at the origin for the constituent Gaussian fields.
  • lkc_ambient (np.array or None): The values for the Lipschitz–Killing Curvatures of the ambient space.
Chi2(dim, dof=2, sigma=1.0, mu=1.0, nu=1.0, lkc_ambient=None)
60    def __init__(self, dim, dof=2, sigma=1., mu=1., nu=1., lkc_ambient=None):
61        super().__init__(dim, name='Isotropic Chi²', sigma=sigma, mu=mu, nu=nu, lkc_ambient=lkc_ambient)
62        self.dof = dof
def LKC(self, j, us):
64    def LKC(self, j, us):
65        """Compute the expected values of the Lipschitz–Killing Curvatures of the excursion sets at thresholds `us`, $\mathbb{L}_j(A_u(f))$.
66
67        Parameters
68        ----------
69        j : int
70            Index of the LKC to compute, `0 < j < dim`.
71            
72        us : np.array
73            The thresholds considered for the computation.
74
75        Returns
76        -------
77        lkc : np.array
78            Expected value of LKC evaluated at the thresholds.
79
80        """
81        us /= self.sigma
82        return LKC_Chi2(j, us, self.mu, dof=self.dof, dim=self.dim, lkc_ambient=self.lkc_ambient)/self.lkc_ambient[-1]

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.
def MF(self, j, us):
 84    def MF(self, j, us):
 85        """Compute the expected values of the Minkowski Functionals of the excursion sets at thresholds `us`, $V_j(A_u(f))$.
 86
 87        Parameters
 88        ----------
 89        j : int
 90            Index of the MF to compute, `0 < j < dim`.
 91            
 92        us : np.array
 93            The thresholds considered for the computation.
 94
 95        Returns
 96        -------
 97        mf : np.array
 98            Expected value of MFs evaluated at the thresholds.
 99
100        """
101        us /= self.sigma
102        return 1./comb(self.dim,j)*LKC_Chi2(self.dim-j, us, self.mu, dof=self.dof, dim=self.dim, lkc_ambient=self.lkc_ambient)/self.lkc_ambient[-1]

Compute the expected values of the Minkowski Functionals of the excursion sets at thresholds us, $V_j(A_u(f))$.

Parameters
  • j (int): Index of the MF to compute, 0 < j < dim.
  • us (np.array): The thresholds considered for the computation.
Returns
  • mf (np.array): Expected value of MFs evaluated at the thresholds.
def V0(self, us):
104    def V0(self, us):
105        """Compute the expected values of the first Minkowski Functionals of the excursion sets at thresholds `us`, $V_0(A_u(f))$.
106
107        Parameters
108        ----------
109        us : np.array
110            The thresholds considered for the computation.
111
112        Returns
113        -------
114        v0 : np.array
115            Expected value of $V_0$ evaluated at the thresholds.
116
117        """
118        us /= self.sigma
119        return self.MF(0, us)

Compute the expected values of the first Minkowski Functionals of the excursion sets at thresholds us, $V_0(A_u(f))$.

Parameters
  • us (np.array): The thresholds considered for the computation.
Returns
  • v0 (np.array): Expected value of $V_0$ evaluated at the thresholds.
def V1(self, us):
121    def V1(self, us):
122        """Compute the expected values of the first Minkowski Functionals of the excursion sets at thresholds `us`, $V_1(A_u(f))$.
123
124        Parameters
125        ----------
126        us : np.array
127            The thresholds considered for the computation.
128
129        Returns
130        -------
131        v1 : np.array
132            Expected value of $V_1$ evaluated at the thresholds.
133
134        """
135        us /= self.sigma
136        return self.MF(1, us)

Compute the expected values of the first Minkowski Functionals of the excursion sets at thresholds us, $V_1(A_u(f))$.

Parameters
  • us (np.array): The thresholds considered for the computation.
Returns
  • v1 (np.array): Expected value of $V_1$ evaluated at the thresholds.
def V2(self, us):
138    def V2(self, us):
139        """Compute the expected values of the first Minkowski Functionals of the excursion sets at thresholds `us`, $V_2(A_u(f))$. Only for fields with `dim ≥ 2`.
140
141        Parameters
142        ----------
143        us : np.array
144            The thresholds considered for the computation.
145
146        Returns
147        -------
148        v2 : np.array
149            Expected value of $V_2$ evaluated at the thresholds.
150
151        """
152        assert self.dim>=2, 'V2 is defined only for fields with dim≥2'
153        us /= self.sigma
154        return self.MF(2, us)

Compute the expected values of the first Minkowski Functionals of the excursion sets at thresholds us, $V_2(A_u(f))$. Only for fields with dim ≥ 2.

Parameters
  • us (np.array): The thresholds considered for the computation.
Returns
  • v2 (np.array): Expected value of $V_2$ evaluated at the thresholds.
def V3(self, us):
156    def V3(self, us):
157        """Compute the expected values of the first Minkowski Functionals of the excursion sets at thresholds `us`, $V_3(A_u(f))$. Only for fields with `dim ≥ 3`.
158
159        Parameters
160        ----------
161        us : np.array
162            The thresholds considered for the computation.
163
164        Returns
165        -------
166        v3 : np.array
167            Expected value of $V_3$ evaluated at the thresholds.
168
169        """
170        assert self.dim>=3, 'V3 is defined only for fields with dim≥3'
171        us /= self.sigma
172        return self.MF(3, us)

Compute the expected values of the first Minkowski Functionals of the excursion sets at thresholds us, $V_3(A_u(f))$. Only for fields with dim ≥ 3.

Parameters
  • us (np.array): The thresholds considered for the computation.
Returns
  • v3 (np.array): Expected value of $V_3$ evaluated at the thresholds.