Skip to content

wave

Polarization ¤

Bases: Module

Source code in jimgw/wave.py
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
class Polarization(eqx.Module):

    name: str
    """Object defining a given polarization mode, with utilities to produce
    corresponding tensor in an Earth centric frame.

    Arguments
    ---------
    name : str
        one of 'p' (plus), 'c' (cross), 'x' (vector x), 'y' (vector y), 'b'
        (breathing), or 'l' (longitudinal).
    """
    def __init__(self, name):
        self.name = name.lower()
        if self.name not in KNOWN_POLS:
            e = f"unknown mode '{self.name}'; must be one of: {KNOWN_POLS}"
            raise ValueError(e)

    def tensor_from_basis(self, x: Array, y: Array) -> Array:
        """Constructor to obtain polarization tensor from waveframe basis
        defined by orthonormal vectors (x, y) in arbitrary Cartesian
        coordinates.
        """
        if self.name == 'p':
            return jnp.einsum('i,j->ij', x, x) - jnp.einsum('i,j->ij', y, y)
        elif self.name == 'c':
            return jnp.einsum('i,j->ij', x, y) + jnp.einsum('i,j->ij', y, x)
        elif self.name == 'x':
            z = jnp.cross(x, y)
            return jnp.einsum('i,j->ij', x, z) + jnp.einsum('i,j->ij', z, x)
        elif self.name == 'y':
            z = jnp.cross(x, y)
            return jnp.einsum('i,j->ij', y, z) + jnp.einsum('i,j->ij', z, y)
        elif self.name == 'b':
            return jnp.einsum('i,j->ij', x, x) + jnp.einsum('i,j->ij', y, y)
        elif self.name == 'l':
            z = jnp.cross(x, y)
            return jnp.einsum('i,j->ij', z, z)
        else:
            raise ValueError(f"unrecognized polarization {self.name}")

    def tensor_from_sky(self, ra: float, dec: float, psi: float, gmst: float) -> Array:
        """Computes {name} polarization tensor in celestial
        coordinates from sky location and orientation parameters.

        Arguments
        ---------
        ra : float
            right ascension in radians.
        dec : float
            declination in radians.
        psi : float
            polarization angle in radians.
        gmst : float
            Greenwhich mean standard time (GMST) in radians.

        Returns
        -------
        tensor : array
            3x3 polarization tensor.
        """
        gmst = jnp.mod(gmst, 2*jnp.pi)
        phi = ra - gmst
        theta = jnp.pi / 2 - dec

        u = jnp.array([jnp.cos(phi) * jnp.cos(theta),
                        jnp.cos(theta) * jnp.sin(phi),
                        -jnp.sin(theta)])
        v = jnp.array([-jnp.sin(phi), jnp.cos(phi), 0])
        m = -u * jnp.sin(psi) - v * jnp.cos(psi)
        n = -u * jnp.cos(psi) + v * jnp.sin(psi)

        return self.tensor_from_basis(m, n)

name: str = name.lower() instance-attribute ¤

Object defining a given polarization mode, with utilities to produce corresponding tensor in an Earth centric frame.

Arguments¤

name : str one of 'p' (plus), 'c' (cross), 'x' (vector x), 'y' (vector y), 'b' (breathing), or 'l' (longitudinal).

tensor_from_basis(x, y) ¤

Constructor to obtain polarization tensor from waveframe basis defined by orthonormal vectors (x, y) in arbitrary Cartesian coordinates.

Source code in jimgw/wave.py
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
def tensor_from_basis(self, x: Array, y: Array) -> Array:
    """Constructor to obtain polarization tensor from waveframe basis
    defined by orthonormal vectors (x, y) in arbitrary Cartesian
    coordinates.
    """
    if self.name == 'p':
        return jnp.einsum('i,j->ij', x, x) - jnp.einsum('i,j->ij', y, y)
    elif self.name == 'c':
        return jnp.einsum('i,j->ij', x, y) + jnp.einsum('i,j->ij', y, x)
    elif self.name == 'x':
        z = jnp.cross(x, y)
        return jnp.einsum('i,j->ij', x, z) + jnp.einsum('i,j->ij', z, x)
    elif self.name == 'y':
        z = jnp.cross(x, y)
        return jnp.einsum('i,j->ij', y, z) + jnp.einsum('i,j->ij', z, y)
    elif self.name == 'b':
        return jnp.einsum('i,j->ij', x, x) + jnp.einsum('i,j->ij', y, y)
    elif self.name == 'l':
        z = jnp.cross(x, y)
        return jnp.einsum('i,j->ij', z, z)
    else:
        raise ValueError(f"unrecognized polarization {self.name}")

tensor_from_sky(ra, dec, psi, gmst) ¤

Computes {name} polarization tensor in celestial coordinates from sky location and orientation parameters.

Arguments¤

ra : float right ascension in radians. dec : float declination in radians. psi : float polarization angle in radians. gmst : float Greenwhich mean standard time (GMST) in radians.

Returns¤

tensor : array 3x3 polarization tensor.

Source code in jimgw/wave.py
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
def tensor_from_sky(self, ra: float, dec: float, psi: float, gmst: float) -> Array:
    """Computes {name} polarization tensor in celestial
    coordinates from sky location and orientation parameters.

    Arguments
    ---------
    ra : float
        right ascension in radians.
    dec : float
        declination in radians.
    psi : float
        polarization angle in radians.
    gmst : float
        Greenwhich mean standard time (GMST) in radians.

    Returns
    -------
    tensor : array
        3x3 polarization tensor.
    """
    gmst = jnp.mod(gmst, 2*jnp.pi)
    phi = ra - gmst
    theta = jnp.pi / 2 - dec

    u = jnp.array([jnp.cos(phi) * jnp.cos(theta),
                    jnp.cos(theta) * jnp.sin(phi),
                    -jnp.sin(theta)])
    v = jnp.array([-jnp.sin(phi), jnp.cos(phi), 0])
    m = -u * jnp.sin(psi) - v * jnp.cos(psi)
    n = -u * jnp.cos(psi) + v * jnp.sin(psi)

    return self.tensor_from_basis(m, n)