Skip to content

utils

A set of utility functions and constants used for unit conversions and cosmology as well as some generically useful math functions.

K_CMB2K_RJ(freq)

Convert from K_CMB to K_RJ.

Parameters:

Name Type Description Default
freq float

The observing frequency in Hz.

required

Returns:

Name Type Description
K_CMB2K_RJ float

Conversion factor from K_CMB to K_RJ.

Source code in witch/utils.py
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
@partial(jax.jit, static_argnums=(0,))
def K_CMB2K_RJ(freq: float) -> float:
    """
    Convert from K_CMB to K_RJ.

    Parameters
    ----------
    freq : float
        The observing frequency in Hz.

    Returns
    -------
    K_CMB2K_RJ : float
        Conversion factor from K_CMB to K_RJ.
    """
    x = freq * h / kb / Tcmb
    return jnp.exp(x) * x * x / jnp.expm1(x) ** 2

beam_double_gauss(dr, fwhm1, amp1, fwhm2, amp2)

Helper function to generate a double gaussian beam.

Parameters:

Name Type Description Default
dr float

Pixel size.

required
fwhm1 float

Full width half max of the primary gaussian in the same units as dr.

required
amp1 float

Amplitude of the primary gaussian.

required
fwhm2 float

Full width half max of the secondairy gaussian in the same units as dr.

required
amp2 float

Amplitude of the secondairy gaussian.

required

Returns:

Type Description
beam: Double gaussian beam.
Source code in witch/utils.py
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
def beam_double_gauss(
    dr: float, fwhm1: float, amp1: float, fwhm2: float, amp2: float
) -> jax.Array:
    """
    Helper function to generate a double gaussian beam.

    Parameters
    ----------
    dr : float
        Pixel size.
    fwhm1 : float
        Full width half max of the primary gaussian in the same units as `dr`.
    amp1 : float
        Amplitude of the primary gaussian.
    fwhm2 : float
        Full width half max of the secondairy gaussian in the same units as `dr`.
    amp2 : float
        Amplitude of the secondairy gaussian.

    Returns
    -------
        beam: Double gaussian beam.
    """
    x = jnp.arange(-1.5 * fwhm1 // (dr), 1.5 * fwhm1 // (dr)) * (dr)
    beam_xx, beam_yy = jnp.meshgrid(x, x)
    beam_rr = jnp.sqrt(beam_xx**2 + beam_yy**2)
    beam = amp1 * jnp.exp(-4 * jnp.log(2) * beam_rr**2 / fwhm1**2) + amp2 * jnp.exp(
        -4 * jnp.log(2) * beam_rr**2 / fwhm2**2
    )
    return beam / jnp.sum(beam)

bilinear_interp(x, y, xp, yp, fp)

JAX implementation of bilinear interpolation. Out of bounds values are set to 0. Using the repeated linear interpolation method here, see https://en.wikipedia.org/wiki/Bilinear_interpolation#Repeated_linear_interpolation.

Parameters:

Name Type Description Default
x Array

X values to return interpolated values at.

required
y Array

Y values to return interpolated values at.

required
xp Array

X values to interpolate with, should be 1D. Assumed to be sorted.

required
yp Array

Y values to interpolate with, should be 1D. Assumed to be sorted.

required
fp Array

Functon values at (xp, yp), should have shape (len(xp), len(yp)). Note that if you are using meshgrid, we assume 'ij' indexing.

required

Returns:

Name Type Description
f Array

The interpolated values.

Source code in witch/utils.py
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
@jax.jit
def bilinear_interp(
    x: jax.Array, y: jax.Array, xp: jax.Array, yp: jax.Array, fp: jax.Array
) -> jax.Array:
    """
    JAX implementation of bilinear interpolation.
    Out of bounds values are set to 0.
    Using the repeated linear interpolation method here,
    see https://en.wikipedia.org/wiki/Bilinear_interpolation#Repeated_linear_interpolation.

    Parameters
    ----------
    x : jax.Array
        X values to return interpolated values at.
    y : jax.Array
        Y values to return interpolated values at.
    xp : jax.Array
        X values to interpolate with, should be 1D.
        Assumed to be sorted.
    yp : jax.Array
        Y values to interpolate with, should be 1D.
        Assumed to be sorted.
    fp : jax.Array
        Functon values at `(xp, yp)`, should have shape `(len(xp), len(yp))`.
        Note that if you are using meshgrid, we assume `'ij'` indexing.

    Returns
    -------
    f : jax.Array
        The interpolated values.
    """
    if len(xp.shape) != 1:
        raise ValueError("xp must be 1D")
    if len(yp.shape) != 1:
        raise ValueError("yp must be 1D")
    if fp.shape != xp.shape + yp.shape:
        raise ValueError(
            "Incompatible shapes for fp, xp, yp: %s, %s, %s",
            fp.shape,
            xp.shape,
            yp.shape,
        )

    # Figure out bounds and mapping
    # This breaks if xp, yp is not sorted
    ix = jnp.clip(jnp.searchsorted(xp, x, side="right"), 1, len(xp) - 1)
    iy = jnp.clip(jnp.searchsorted(yp, y, side="right"), 1, len(yp) - 1)
    q_11 = fp[ix - 1, iy - 1]
    q_21 = fp[ix, iy - 1]
    q_12 = fp[ix - 1, iy]
    q_22 = fp[ix, iy]

    # Interpolate in x to start
    denom_x = xp[ix] - xp[ix - 1]
    dx_1 = x - xp[ix - 1]
    dx_2 = xp[ix] - x
    f_xy1 = (dx_2 * q_11 + dx_1 * q_21) / denom_x
    f_xy2 = (dx_2 * q_12 + dx_1 * q_22) / denom_x

    # Now do y as well
    denom_y = yp[iy] - yp[iy - 1]
    dy_1 = y - yp[iy - 1]
    dy_2 = yp[iy] - y
    f = (dy_2 * f_xy1 + dy_1 * f_xy2) / denom_y

    # Zero out the out of bounds values
    f = jnp.where((x < xp[0]) + (x > xp[-1]) + (y < yp[0]) + (y > yp[-1]), 0.0, f)

    return f

bin_map(data, pixsize)

Bins data radially.

Parameters:

Name Type Description Default
data ArrayLike

Data to be radially binned

required
pixsize float

Pixel spacing for data

required

Returns:

Name Type Description
rs array

Left bin edges

bin1d array

Mean of pixels in bin

var1d array

Variance of pixels in bin

Source code in witch/utils.py
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
def bin_map(data: ArrayLike, pixsize: float) -> tuple[np.array, np.array, np.array]:
    """
    Bins data radially.

    Parameters
    ----------
    data : ArrayLike
        Data to be radially binned
    pixsize : float
        Pixel spacing for data

    Returns
    -------
    rs : np.array
        Left bin edges
    bin1d : np.array
        Mean of pixels in bin
    var1d : np.array
        Variance of pixels in bin
    """
    x = np.linspace(
        -data.shape[1] / 2 * pixsize, data.shape[1] / 2 * pixsize, data.shape[1]
    )
    y = np.linspace(
        -data.shape[0] / 2 * pixsize, data.shape[0] / 2 * pixsize, data.shape[0]
    )

    X, Y = np.meshgrid(x, y)
    R = np.sqrt(X**2 + Y**2)  # TODO: miscentering?

    rs = np.arange(0, np.amax(R), pixsize)
    rs = np.append(rs, 999999)
    bin1d = np.zeros(len(rs) - 1)
    var1d = np.zeros(len(rs) - 1)

    for k in range(len(rs) - 1):
        pixels = [
            data[i, j]
            for i in range(len(y))
            for j in range(len(x))
            if rs[k] < R[i, j] <= rs[k + 1]
        ]
        if len(pixels) == 0:
            bin1d[k] = 0
            var1d[k] = 0
        else:
            bin1d[k] = np.mean(pixels)
            var1d[k] = np.var(pixels)
    rs = rs[:-1]

    return rs, bin1d, var1d

fft_conv(image, kernel)

Perform a convolution using FFTs for speed with jax.

Parameters:

Name Type Description Default
image ArrayLike

Data to be convolved.

required
kernel ArrayLike

Convolution kernel.

required

Returns:

Name Type Description
convolved_map Array

Image convolved with kernel.

Source code in witch/utils.py
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
@jax.jit
def fft_conv(image: ArrayLike, kernel: ArrayLike) -> jax.Array:
    """
    Perform a convolution using FFTs for speed with jax.

    Parameters
    ----------
    image : ArrayLike
        Data to be convolved.
    kernel : ArrayLike
        Convolution kernel.

    Returns
    -------
    convolved_map : jax.Array
        Image convolved with kernel.
    """
    Fmap = jnp.fft.fft2(jnp.fft.fftshift(image))
    Fkernel = jnp.fft.fft2(jnp.fft.fftshift(kernel))
    convolved_map = jnp.fft.fftshift(jnp.real(jnp.fft.ifft2(Fmap * Fkernel)))

    return convolved_map

fft_deconv(image, kernel)

Perform a convolution using FFTs for speed with jax.

Parameters:

Name Type Description Default
image ArrayLike

Data to be convolved.

required
kernel ArrayLike

Convolution kernel.

required

Returns:

Name Type Description
convolved_map Array

Image convolved with kernel.

Source code in witch/utils.py
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
@jax.jit
def fft_deconv(image: ArrayLike, kernel: ArrayLike) -> jax.Array:
    """
    Perform a convolution using FFTs for speed with jax.

    Parameters
    ----------
    image : ArrayLike
        Data to be convolved.
    kernel : ArrayLike
        Convolution kernel.

    Returns
    -------
    convolved_map : jax.Array
        Image convolved with kernel.
    """
    Fmap = jnp.fft.fft2(jnp.fft.fftshift(image))
    Fkernel = jnp.fft.fft2(jnp.fft.fftshift(kernel))
    convolved_map = jnp.fft.fftshift(jnp.real(jnp.fft.ifft2(Fmap / Fkernel)))

    return convolved_map

get_da(z)

Get factor to convert from arcseconds to MPc.

Parameters:

Name Type Description Default
z ArrayLike

The redshift at which to compute the factor.

required

Returns:

Name Type Description
da Array

Conversion factor from arcseconds to MPc.

Source code in witch/utils.py
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
@jax.jit
def get_da(z: ArrayLike) -> jax.Array:
    """
    Get factor to convert from arcseconds to MPc.

    Parameters
    ----------
    z : ArrayLike
        The redshift at which to compute the factor.

    Returns
    -------
    da : jax.Array
        Conversion factor from arcseconds to MPc.
    """
    return jnp.interp(z, zline, daline)

get_hz(z)

Get the dimensionless hubble constant, h, at a given redshift.

Parameters:

Name Type Description Default
z ArrayLike

The redshift at which to compute h.

required

Returns:

Name Type Description
hz Array

h at the given z.

Source code in witch/utils.py
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
@jax.jit
def get_hz(z: ArrayLike) -> jax.Array:
    """
    Get the dimensionless hubble constant, h, at a given redshift.

    Parameters
    ----------
    z : ArrayLike
        The redshift at which to compute h.

    Returns
    -------
    hz : jax.Array
        h at the given z.
    """
    return jnp.interp(z, zline, hzline)

get_nz(z)

Get the critical density at a given redshift.

Parameters:

Name Type Description Default
z ArrayLike

The redshift at which to compute the critical density.

required

Returns:

Name Type Description
nz Array

Critical density at the given z. This is in units of solar masses per cubic Mpc.

Source code in witch/utils.py
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
@jax.jit
def get_nz(z: ArrayLike) -> jax.Array:
    """
    Get the critical density at a given redshift.

    Parameters
    ----------
    z : ArrayLike
        The redshift at which to compute the critical density.

    Returns
    -------
    nz : jax.Array
        Critical density at the given z.
        This is in units of solar masses per cubic Mpc.
    """
    return jnp.interp(z, zline, nzline)

tod_hi_pass(tod, N_filt)

High pass a tod with a tophat

Parameters:

Name Type Description Default
tod Array

TOD to high pass.

required
N_filt int

N_filt of tophat.

required

Returns:

Name Type Description
tod_filtered Array

High pass filtered TOD

Source code in witch/utils.py
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
@partial(jax.jit, static_argnums=(1,))
def tod_hi_pass(tod: jax.Array, N_filt: int) -> jax.Array:
    """
    High pass a tod with a tophat

    Parameters
    ----------
    tod : jax.Array
        TOD to high pass.
    N_filt : int
        N_filt of tophat.

    Returns
    -------
    tod_filtered : jax.Array
        High pass filtered TOD
    """
    mask = jnp.ones(tod.shape)
    mask = mask.at[..., :N_filt].set(0.0)

    ## apply the filter in fourier space
    Ftod = jnp.fft.fft(tod)
    Ftod_filtered = Ftod * mask
    tod_filtered = jnp.fft.ifft(Ftod_filtered).real
    return tod_filtered

y2K_CMB(freq, Te)

Convert from compton y to K_CMB.

Parameters:

Name Type Description Default
freq float

The observing frequency in Hz.

required
Te float

Electron temperature

required

Returns:

Name Type Description
y2K_CMB float

Conversion factor from compton y to K_CMB.

Source code in witch/utils.py
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
@partial(jax.jit, static_argnums=(0, 1))
def y2K_CMB(freq: float, Te: float) -> float:
    """
    Convert from compton y to K_CMB.

    Parameters
    ----------
    freq : float
        The observing frequency in Hz.
    Te : float
        Electron temperature

    Returns
    -------
    y2K_CMB : float
        Conversion factor from compton y to K_CMB.
    """
    x = freq * h / kb / Tcmb
    xt = x / jnp.tanh(0.5 * x)
    st = x / jnp.sinh(0.5 * x)
    # fmt:off
    Y0 = -4.0 + xt
    Y1 = (-10.0
        + ((47.0 / 2.0) + (-(42.0 / 5.0) + (7.0 / 10.0) * xt) * xt) * xt
        + st * st * (-(21.0 / 5.0) + (7.0 / 5.0) * xt)
    )
    Y2 = ((-15.0 / 2.0)
        + ((1023.0 / 8.0) + ((-868.0 / 5.0) + ((329.0 / 5.0) + ((-44.0 / 5.0) + (11.0 / 30.0) * xt) * xt) * xt) * xt) * xt
        + ((-434.0 / 5.0) + ((658.0 / 5.0) + ((-242.0 / 5.0) + (143.0 / 30.0) * xt) * xt) * xt
        + (-(44.0 / 5.0) + (187.0 / 60.0) * xt) * (st * st)) * st * st
    )
    Y3 = ((15.0 / 2.0)
        + ((2505.0 / 8.0) + ((-7098.0 / 5.0) + ((14253.0 / 10.0) + ((-18594.0 / 35.0) 
         + ((12059.0 / 140.0) + ((-128.0 / 21.0) + (16.0 / 105.0) * xt) * xt) * xt) * xt) * xt) * xt) * xt
        + (((-7098.0 / 10.0) + ((14253.0 / 5.0) + ((-102267.0 / 35.0) + ((156767.0 / 140.0)
         + ((-1216.0 / 7.0) + (64.0 / 7.0) * xt) * xt) * xt) * xt) * xt)
         + (((-18594.0 / 35.0) + ((205003.0 / 280.0) + ((-1920.0 / 7.0) + (1024.0 / 35.0) * xt) * xt) * xt)
          + ((-544.0 / 21.0) + (992.0 / 105.0) * xt) * st * st) * st * st) * st * st
    )
    Y4 = ((-135.0 / 32.0)
        + ((30375.0 / 128.0) + ((-62391.0 / 10.0) + ((614727.0 / 40.0) + ((-124389.0 / 10.0) + ((355703.0 / 80.0) + ((-16568.0 / 21.0)
         + ((7516.0 / 105.0) + ((-22.0 / 7.0) + (11.0 / 210.0) * xt) * xt) * xt) * xt) * xt) * xt) * xt) * xt) * xt
        + ((-62391.0 / 20.0) + ((614727.0 / 20.0) + ((-1368279.0 / 20.0) + ((4624139.0 / 80.0) + ((-157396.0 / 7.0) + ((30064.0 / 7.0)
         + ((-2717.0 / 7.0) + (2761.0 / 210.0) * xt) * xt) * xt) * xt) * xt) * xt) * xt
         + ((-124389.0 / 10.0)
          + ((6046951.0 / 160.0) + ((-248520.0 / 7.0) + ((481024.0 / 35.0) + ((-15972.0 / 7.0) + (18689.0 / 140.0) * xt) * xt) * xt) * xt) * xt
          + ((-70414.0 / 21.0) + ((465992.0 / 105.0) + ((-11792.0 / 7.0) + (19778.0 / 105.0) * xt) * xt) * xt
           + ((-682.0 / 7.0) + (7601.0 / 210.0) * xt) * st * st) * st * st) * st * st) * st * st
    )
    # fmt:on
    factor = Y0 + (Te / me) * (
        Y1 + (Te / me) * (Y2 + (Te / me) * (Y3 + (Te / me) * Y4))
    )
    return factor * Tcmb

y2K_RJ(freq, Te)

Convert from compton y to K_RJ.

Parameters:

Name Type Description Default
freq float

The observing frequency in Hz.

required
Te float

Electron temperature

required

Returns:

Name Type Description
y2K_RJ float

Conversion factor from compton y to K_RJ.

Source code in witch/utils.py
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
@partial(jax.jit, static_argnums=(0, 1))
def y2K_RJ(freq: float, Te: float) -> float:
    """
    Convert from compton y to K_RJ.

    Parameters
    ----------
    freq : float
        The observing frequency in Hz.
    Te : float
        Electron temperature

    Returns
    -------
    y2K_RJ : float
        Conversion factor from compton y to K_RJ.
    """
    factor = y2K_CMB(freq, Te)
    return factor * K_CMB2K_RJ(freq)

y2uK_CMB(freq, Te)

Convert from compton y to uK_CMB.

Parameters:

Name Type Description Default
freq float

The observing frequency in Hz.

required
Te float

Electron temperature

required

Returns:

Name Type Description
y2uK_RJ float

Conversion factor from compton y to uK_CMB.

Source code in witch/utils.py
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
@partial(jax.jit, static_argnums=(0, 1))
def y2uK_CMB(freq: float, Te: float) -> float:
    """
    Convert from compton y to uK_CMB.

    Parameters
    ----------
    freq : float
        The observing frequency in Hz.
    Te : float
        Electron temperature

    Returns
    -------
    y2uK_RJ : float
        Conversion factor from compton y to uK_CMB.
    """
    factor = y2K_CMB(freq, Te)
    return 1e6 * factor