Skip to content

Instrument Model

For each Observation (see Data Containers), a unique instrument model is specified. This includes the set of detectors, their properties, and other metadata about the overall telescope.

toast.instrument.Site

Bases: object

Site base class.

Parameters:

Name Type Description Default
name str

Site name

required
uid int

Unique identifier. If not specified, constructed from a hash of the site name.

None
Source code in toast/instrument.py
 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
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
class Site(object):
    """Site base class.

    Args:
        name (str):  Site name
        uid (int):  Unique identifier.  If not specified, constructed from a hash
            of the site name.

    """

    def __init__(self, name, uid=None):
        self.name = name
        self.uid = uid
        if self.uid is None:
            self.uid = name_UID(self.name)

    def _position(self, times):
        raise NotImplementedError("Derived class must implement _position()")

    def position(self, times):
        """Get the site position in solar system barycentric cartesian vectors.

        Given timestamps in POSIX seconds since 1970 (UTC), return the position as
        solar system coordinates.

        Args:
            times (array):  The timestamps.

        Returns:
            (array):  The position vectors.

        """
        return self._position(times)

    def _velocity(self, times):
        raise NotImplementedError("Derived class must implement _velocity()")

    def velocity(self, times):
        """Get the site velocity in solar system barycentric cartesian vectors.

        Given timestamps in POSIX seconds since 1970 (UTC), return the velocity as
        quaternions in solar system barycentric coordinates.

        Args:
            times (array):  The timestamps.

        Returns:
            (array):  The velocity vectors.

        """
        return self._velocity(times)

    def position_velocity(self, times):
        """Get the site position and velocity.

        Convenience function to simultaneously return the position and velocity.

        Args:
            times (array):  The timestamps.

        Returns:
            (tuple):  The position and velocity arrays of vectors.

        """
        if hasattr(self, "_position_velocity"):
            return self._position_velocity(times)
        else:
            p = self._position(times)
            v = self._velocity(times)
        return (p, v)

    def __repr__(self):
        value = "<Site '{}' : uid = {}>".format(self.name, self.uid)
        return value

    def __eq__(self, other):
        if self.name != other.name:
            return False
        if self.uid != other.uid:
            return False
        return True

    def __ne__(self, other):
        return not self.__eq__(other)

name = name instance-attribute

uid = uid instance-attribute

__eq__(other)

Source code in toast/instrument.py
117
118
119
120
121
122
def __eq__(self, other):
    if self.name != other.name:
        return False
    if self.uid != other.uid:
        return False
    return True

__init__(name, uid=None)

Source code in toast/instrument.py
52
53
54
55
56
def __init__(self, name, uid=None):
    self.name = name
    self.uid = uid
    if self.uid is None:
        self.uid = name_UID(self.name)

__ne__(other)

Source code in toast/instrument.py
124
125
def __ne__(self, other):
    return not self.__eq__(other)

__repr__()

Source code in toast/instrument.py
113
114
115
def __repr__(self):
    value = "<Site '{}' : uid = {}>".format(self.name, self.uid)
    return value

_position(times)

Source code in toast/instrument.py
58
59
def _position(self, times):
    raise NotImplementedError("Derived class must implement _position()")

_velocity(times)

Source code in toast/instrument.py
76
77
def _velocity(self, times):
    raise NotImplementedError("Derived class must implement _velocity()")

position(times)

Get the site position in solar system barycentric cartesian vectors.

Given timestamps in POSIX seconds since 1970 (UTC), return the position as solar system coordinates.

Parameters:

Name Type Description Default
times array

The timestamps.

required

Returns:

Type Description
array

The position vectors.

Source code in toast/instrument.py
61
62
63
64
65
66
67
68
69
70
71
72
73
74
def position(self, times):
    """Get the site position in solar system barycentric cartesian vectors.

    Given timestamps in POSIX seconds since 1970 (UTC), return the position as
    solar system coordinates.

    Args:
        times (array):  The timestamps.

    Returns:
        (array):  The position vectors.

    """
    return self._position(times)

position_velocity(times)

Get the site position and velocity.

Convenience function to simultaneously return the position and velocity.

Parameters:

Name Type Description Default
times array

The timestamps.

required

Returns:

Type Description
tuple

The position and velocity arrays of vectors.

Source code in toast/instrument.py
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
def position_velocity(self, times):
    """Get the site position and velocity.

    Convenience function to simultaneously return the position and velocity.

    Args:
        times (array):  The timestamps.

    Returns:
        (tuple):  The position and velocity arrays of vectors.

    """
    if hasattr(self, "_position_velocity"):
        return self._position_velocity(times)
    else:
        p = self._position(times)
        v = self._velocity(times)
    return (p, v)

velocity(times)

Get the site velocity in solar system barycentric cartesian vectors.

Given timestamps in POSIX seconds since 1970 (UTC), return the velocity as quaternions in solar system barycentric coordinates.

Parameters:

Name Type Description Default
times array

The timestamps.

required

Returns:

Type Description
array

The velocity vectors.

Source code in toast/instrument.py
79
80
81
82
83
84
85
86
87
88
89
90
91
92
def velocity(self, times):
    """Get the site velocity in solar system barycentric cartesian vectors.

    Given timestamps in POSIX seconds since 1970 (UTC), return the velocity as
    quaternions in solar system barycentric coordinates.

    Args:
        times (array):  The timestamps.

    Returns:
        (array):  The velocity vectors.

    """
    return self._velocity(times)

toast.instrument.GroundSite

Bases: Site

Site on the Earth.

This represents a fixed location on the Earth.

Parameters:

Name Type Description Default
name str

Site name

required
lat Quantity

Site latitude.

required
lon Quantity

Site longitude.

required
alt Quantity

Site altitude.

required
uid int

Unique identifier. If not specified, constructed from a hash of the site name.

None
weather Weather

Weather information for this site.

None
Source code in toast/instrument.py
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
class GroundSite(Site):
    """Site on the Earth.

    This represents a fixed location on the Earth.

    Args:
        name (str):  Site name
        lat (Quantity):  Site latitude.
        lon (Quantity):  Site longitude.
        alt (Quantity):  Site altitude.
        uid (int):  Unique identifier.  If not specified, constructed from a hash
            of the site name.
        weather (Weather):  Weather information for this site.
    """

    def __init__(self, name, lat, lon, alt, uid=None, weather=None):
        super().__init__(name, uid)
        self.earthloc = coord.EarthLocation.from_geodetic(lon, lat, height=alt)
        self.weather = weather

    def __repr__(self):
        value = "<GroundSite '{}' : uid = {}, lon = {}, lat = {}, alt = {}, weather = {}>".format(
            self.name,
            self.uid,
            self.earthloc.lon,
            self.earthloc.lat,
            self.earthloc.height,
            self.weather,
        )
        return value

    def __eq__(self, other):
        if self.name != other.name:
            return False
        if self.uid != other.uid:
            return False
        if not np.isclose(other.earthloc.lon, self.earthloc.lon):
            return False
        if not np.isclose(other.earthloc.lat, self.earthloc.lat):
            return False
        if not np.isclose(other.earthloc.height, self.earthloc.height):
            return False
        if self.weather != other.weather:
            return False
        return True

    def _position_velocity(self, times):
        # Compute data at 10 second intervals and interpolate.  If the timestamps are
        # more coarsely sampled than that, just compute those times directly.
        sparse_incr = 10.0
        do_interp = True
        if len(times) < 100 or (times[1] - times[0]) > sparse_incr:
            do_interp = False

        if do_interp:
            n_sparse = int((times[-1] - times[0]) / sparse_incr)
            sparse_times = np.linspace(times[0], times[-1], num=n_sparse, endpoint=True)
        else:
            n_sparse = len(times)
            sparse_times = times
        pos_x = np.zeros(n_sparse, np.float64)
        pos_y = np.zeros(n_sparse, np.float64)
        pos_z = np.zeros(n_sparse, np.float64)
        vel_x = np.zeros(n_sparse, np.float64)
        vel_y = np.zeros(n_sparse, np.float64)
        vel_z = np.zeros(n_sparse, np.float64)
        for i, t in enumerate(sparse_times):
            atime = astime.Time(t, format="unix")
            p, v = coord.get_body_barycentric_posvel("earth", atime)
            # FIXME:  apply translation from earth center to earth location.
            # itrs = self.earthloc.get_itrs(obstime)
            pm = p.xyz.to_value(u.kilometer)
            vm = v.xyz.to_value(u.kilometer / u.second)
            pos_x[i] = pm[0]
            pos_y[i] = pm[1]
            pos_z[i] = pm[2]
            vel_x[i] = vm[0]
            vel_y[i] = vm[1]
            vel_z[i] = vm[2]

        if do_interp:
            pos_x = np.interp(times, sparse_times, pos_x)
            pos_y = np.interp(times, sparse_times, pos_y)
            pos_z = np.interp(times, sparse_times, pos_z)
            vel_x = np.interp(times, sparse_times, vel_x)
            vel_y = np.interp(times, sparse_times, vel_y)
            vel_z = np.interp(times, sparse_times, vel_z)
        pos = np.stack([pos_x, pos_y, pos_z], axis=-1)
        vel = np.stack([vel_x, vel_y, vel_z], axis=-1)
        return pos, vel

    def _position(self, times):
        p, v = self._position_velocity(times)
        return p

    def _velocity(self, times):
        p, v = self._position_velocity(times)
        return v

earthloc = coord.EarthLocation.from_geodetic(lon, lat, height=alt) instance-attribute

weather = weather instance-attribute

__eq__(other)

Source code in toast/instrument.py
159
160
161
162
163
164
165
166
167
168
169
170
171
172
def __eq__(self, other):
    if self.name != other.name:
        return False
    if self.uid != other.uid:
        return False
    if not np.isclose(other.earthloc.lon, self.earthloc.lon):
        return False
    if not np.isclose(other.earthloc.lat, self.earthloc.lat):
        return False
    if not np.isclose(other.earthloc.height, self.earthloc.height):
        return False
    if self.weather != other.weather:
        return False
    return True

__init__(name, lat, lon, alt, uid=None, weather=None)

Source code in toast/instrument.py
143
144
145
146
def __init__(self, name, lat, lon, alt, uid=None, weather=None):
    super().__init__(name, uid)
    self.earthloc = coord.EarthLocation.from_geodetic(lon, lat, height=alt)
    self.weather = weather

__repr__()

Source code in toast/instrument.py
148
149
150
151
152
153
154
155
156
157
def __repr__(self):
    value = "<GroundSite '{}' : uid = {}, lon = {}, lat = {}, alt = {}, weather = {}>".format(
        self.name,
        self.uid,
        self.earthloc.lon,
        self.earthloc.lat,
        self.earthloc.height,
        self.weather,
    )
    return value

_position(times)

Source code in toast/instrument.py
219
220
221
def _position(self, times):
    p, v = self._position_velocity(times)
    return p

_position_velocity(times)

Source code in toast/instrument.py
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
def _position_velocity(self, times):
    # Compute data at 10 second intervals and interpolate.  If the timestamps are
    # more coarsely sampled than that, just compute those times directly.
    sparse_incr = 10.0
    do_interp = True
    if len(times) < 100 or (times[1] - times[0]) > sparse_incr:
        do_interp = False

    if do_interp:
        n_sparse = int((times[-1] - times[0]) / sparse_incr)
        sparse_times = np.linspace(times[0], times[-1], num=n_sparse, endpoint=True)
    else:
        n_sparse = len(times)
        sparse_times = times
    pos_x = np.zeros(n_sparse, np.float64)
    pos_y = np.zeros(n_sparse, np.float64)
    pos_z = np.zeros(n_sparse, np.float64)
    vel_x = np.zeros(n_sparse, np.float64)
    vel_y = np.zeros(n_sparse, np.float64)
    vel_z = np.zeros(n_sparse, np.float64)
    for i, t in enumerate(sparse_times):
        atime = astime.Time(t, format="unix")
        p, v = coord.get_body_barycentric_posvel("earth", atime)
        # FIXME:  apply translation from earth center to earth location.
        # itrs = self.earthloc.get_itrs(obstime)
        pm = p.xyz.to_value(u.kilometer)
        vm = v.xyz.to_value(u.kilometer / u.second)
        pos_x[i] = pm[0]
        pos_y[i] = pm[1]
        pos_z[i] = pm[2]
        vel_x[i] = vm[0]
        vel_y[i] = vm[1]
        vel_z[i] = vm[2]

    if do_interp:
        pos_x = np.interp(times, sparse_times, pos_x)
        pos_y = np.interp(times, sparse_times, pos_y)
        pos_z = np.interp(times, sparse_times, pos_z)
        vel_x = np.interp(times, sparse_times, vel_x)
        vel_y = np.interp(times, sparse_times, vel_y)
        vel_z = np.interp(times, sparse_times, vel_z)
    pos = np.stack([pos_x, pos_y, pos_z], axis=-1)
    vel = np.stack([vel_x, vel_y, vel_z], axis=-1)
    return pos, vel

_velocity(times)

Source code in toast/instrument.py
223
224
225
def _velocity(self, times):
    p, v = self._position_velocity(times)
    return v

toast.instrument.SpaceSite

Bases: Site

Site with no atmosphere.

This represents a location beyond the Earth's atmosphere. In practice, this should be sub-classed for real satellite experiments.

Parameters:

Name Type Description Default
name str

Site name

required
uid int

Unique identifier. If not specified, constructed from a hash of the site name.

None
Source code in toast/instrument.py
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
class SpaceSite(Site):
    """Site with no atmosphere.

    This represents a location beyond the Earth's atmosphere.  In practice, this
    should be sub-classed for real satellite experiments.

    Args:
        name (str):  Site name
        uid (int):  Unique identifier.  If not specified, constructed from a hash
            of the site name.

    """

    def __init__(self, name, uid=None):
        super().__init__(name, uid)

    def __repr__(self):
        value = "<SpaceSite '{}' : uid = {}>".format(self.name, self.uid)
        return value

    def _position_velocity(self, times):
        # Compute data at 10 minute intervals and interpolate.  If the timestamps are
        # more coarsely sampled than that, just compute those times directly.
        sparse_incr = 600.0
        do_interp = True
        if len(times) < 100 or (times[1] - times[0]) > sparse_incr:
            do_interp = False

        if do_interp:
            n_sparse = 1 + int((times[-1] - times[0]) / sparse_incr)
            sparse_times = np.linspace(times[0], times[-1], num=n_sparse, endpoint=True)
        else:
            n_sparse = len(times)
            sparse_times = times
        pos_x = np.zeros(n_sparse, np.float64)
        pos_y = np.zeros(n_sparse, np.float64)
        pos_z = np.zeros(n_sparse, np.float64)
        vel_x = np.zeros(n_sparse, np.float64)
        vel_y = np.zeros(n_sparse, np.float64)
        vel_z = np.zeros(n_sparse, np.float64)
        for i, t in enumerate(sparse_times):
            atime = astime.Time(t, format="unix")
            # Get the satellite position and velocity in the equatorial frame (ICRS)
            p, v = coord.get_body_barycentric_posvel("earth", atime)
            # FIXME:  apply translation from earth center to L2.
            pm = p.xyz.to_value(u.kilometer)
            vm = v.xyz.to_value(u.kilometer / u.second)
            pos_x[i] = pm[0]
            pos_y[i] = pm[1]
            pos_z[i] = pm[2]
            vel_x[i] = vm[0]
            vel_y[i] = vm[1]
            vel_z[i] = vm[2]

        if do_interp:
            pos_x = np.interp(times, sparse_times, pos_x)
            pos_y = np.interp(times, sparse_times, pos_y)
            pos_z = np.interp(times, sparse_times, pos_z)
            vel_x = np.interp(times, sparse_times, vel_x)
            vel_y = np.interp(times, sparse_times, vel_y)
            vel_z = np.interp(times, sparse_times, vel_z)
        pos = np.stack([pos_x, pos_y, pos_z], axis=-1)
        vel = np.stack([vel_x, vel_y, vel_z], axis=-1)
        return pos, vel

    def _position(self, times):
        p, v = self._position_velocity(times)
        return p

    def _velocity(self, times):
        p, v = self._position_velocity(times)
        return v

__init__(name, uid=None)

Source code in toast/instrument.py
241
242
def __init__(self, name, uid=None):
    super().__init__(name, uid)

__repr__()

Source code in toast/instrument.py
244
245
246
def __repr__(self):
    value = "<SpaceSite '{}' : uid = {}>".format(self.name, self.uid)
    return value

_position(times)

Source code in toast/instrument.py
293
294
295
def _position(self, times):
    p, v = self._position_velocity(times)
    return p

_position_velocity(times)

Source code in toast/instrument.py
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
def _position_velocity(self, times):
    # Compute data at 10 minute intervals and interpolate.  If the timestamps are
    # more coarsely sampled than that, just compute those times directly.
    sparse_incr = 600.0
    do_interp = True
    if len(times) < 100 or (times[1] - times[0]) > sparse_incr:
        do_interp = False

    if do_interp:
        n_sparse = 1 + int((times[-1] - times[0]) / sparse_incr)
        sparse_times = np.linspace(times[0], times[-1], num=n_sparse, endpoint=True)
    else:
        n_sparse = len(times)
        sparse_times = times
    pos_x = np.zeros(n_sparse, np.float64)
    pos_y = np.zeros(n_sparse, np.float64)
    pos_z = np.zeros(n_sparse, np.float64)
    vel_x = np.zeros(n_sparse, np.float64)
    vel_y = np.zeros(n_sparse, np.float64)
    vel_z = np.zeros(n_sparse, np.float64)
    for i, t in enumerate(sparse_times):
        atime = astime.Time(t, format="unix")
        # Get the satellite position and velocity in the equatorial frame (ICRS)
        p, v = coord.get_body_barycentric_posvel("earth", atime)
        # FIXME:  apply translation from earth center to L2.
        pm = p.xyz.to_value(u.kilometer)
        vm = v.xyz.to_value(u.kilometer / u.second)
        pos_x[i] = pm[0]
        pos_y[i] = pm[1]
        pos_z[i] = pm[2]
        vel_x[i] = vm[0]
        vel_y[i] = vm[1]
        vel_z[i] = vm[2]

    if do_interp:
        pos_x = np.interp(times, sparse_times, pos_x)
        pos_y = np.interp(times, sparse_times, pos_y)
        pos_z = np.interp(times, sparse_times, pos_z)
        vel_x = np.interp(times, sparse_times, vel_x)
        vel_y = np.interp(times, sparse_times, vel_y)
        vel_z = np.interp(times, sparse_times, vel_z)
    pos = np.stack([pos_x, pos_y, pos_z], axis=-1)
    vel = np.stack([vel_x, vel_y, vel_z], axis=-1)
    return pos, vel

_velocity(times)

Source code in toast/instrument.py
297
298
299
def _velocity(self, times):
    p, v = self._position_velocity(times)
    return v

toast.instrument.Bandpass

Bases: object

Class that contains the bandpass information for an entire focalplane.

Source code in toast/instrument.py
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
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
390
391
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
422
423
424
425
426
427
428
429
430
431
432
433
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
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
class Bandpass(object):
    """Class that contains the bandpass information for an entire focalplane."""

    @function_timer
    def __init__(self, bandcenters, bandwidths, nstep=101):
        """All units in GHz

        Args :
            bandcenters(dict) : Dictionary of bandpass centers
            bandwidths(dict) : Dictionary of bandpass widths
            nstep(int) : Number of interplation steps to use in `convolve()`
        """
        self.nstep = nstep
        self.dets = []
        self._fmin = {}
        self._fmax = {}
        for name in bandcenters:
            self.dets.append(name)
            center = bandcenters[name]
            width = bandwidths[name]
            self._fmin[name] = center - 0.5 * width
            self._fmax[name] = center + 0.5 * width
        # The interpolated bandpasses will be cached as needed
        self._fmin_tot = None
        self._fmax_tot = None
        self._freqs = {}
        self._bandpass = {}
        self._kcmb2jysr = {}
        self._kcmb2krj = {}
        self._kcmb2w = {}

    @function_timer
    def get_range(self, det=None):
        """Return the maximum range of frequencies needed for convolution."""
        if det is not None:
            return self._fmin[det], self._fmax[det]
        elif self._fmin_tot is None:
            self._fmin_tot = min(self._fmin.values())
            self._fmax_tot = max(self._fmax.values())
        return self._fmin_tot, self._fmax_tot

    @function_timer
    def center_frequency(self, det, alpha=-1):
        """Return the effective central frequency for a given spectral index"""

        # Which delta function bandpass would produce the same flux density
        freqs = self.freqs(det)
        if alpha == 0:
            # The equation is singular at alpha == 0. Evaluate it on both sides
            # and return the average
            delta = 1e-6
            alpha1 = alpha - delta
            eff1 = self.convolve(det, freqs, freqs.to_value(u.Hz) ** alpha1) ** (
                1 / alpha1
            )
            alpha2 = alpha + delta
            eff2 = self.convolve(det, freqs, freqs.to_value(u.Hz) ** alpha2) ** (
                1 / alpha2
            )
            eff = 0.5 * (eff1 + eff2)
        else:
            # Very simple closed form
            eff = self.convolve(det, freqs, freqs.to_value(u.Hz) ** alpha) ** (
                1 / alpha
            )

        return eff * u.Hz

    @function_timer
    def optical_loading(self, det, T):
        """Return the optical loading of a blackbody source.
        Assumes a diffraction-limited, single-moded polarimeter
        and perfect optical efficiency
        arXiv:1806.04316

        Args:
            det(str) : detector name
            T(float) : Source temperature in Kelvin

        Returns:
            (float) : The optical loading in Watts
        """

        bandpass = self.bandpass(det)
        # Normalize the bandpass to peak at 1
        bandpass = bandpass / np.amax(bandpass)
        freqs = self.freqs(det).to_value(u.Hz)

        # Power spectral density
        S = h * freqs / (np.exp(h * freqs / k / T) - 1)

        # Integrate over frequency to get power
        power = integrate_simpson(freqs, S * bandpass)

        return power

    @function_timer
    def _get_unit_conversion_coefficients(self, det):
        """Compute and cache the unit conversion coefficients for one detector"""

        if (
            det not in self._kcmb2jysr
            or det not in self._kcmb2krj
            or det not in self._kcmb2w
        ):
            # The calculation is a copy from the Hildebrandt and Macias-Perez IDL module for Planck

            nu_cmb = k * TCMB / h
            alpha = 2 * k**3 * TCMB**2 / h**2 / c**2

            cfreq = self.center_frequency(det).to_value(u.Hz)
            freqs = self.freqs(det).to_value(u.Hz)
            bandpass = self.bandpass(det)

            x = freqs / nu_cmb
            db_dt = alpha * x**4 * np.exp(x) / (np.exp(x) - 1) ** 2
            db_dt_rj = 2 * freqs**2 * k / c**2

            self._kcmb2jysr[det] = (
                1e26
                * integrate_simpson(freqs, db_dt * bandpass)
                / integrate_simpson(freqs, cfreq / freqs * bandpass)
            )
            self._kcmb2krj[det] = integrate_simpson(
                freqs, db_dt * bandpass
            ) / integrate_simpson(freqs, db_dt_rj * bandpass)

            # K_CMB->W conversion is from the BoloCalc paper, arXiv:1806.04316
            bandpass = bandpass / np.amax(bandpass)
            self._kcmb2w[det] = integrate_simpson(
                freqs,
                k * (x / (np.exp(x) - 1)) ** 2 * np.exp(x) * bandpass,
            )

        return

    @function_timer
    def freqs(self, det):
        if det not in self._freqs:
            fmin = self._fmin[det].to_value(u.Hz)
            fmax = self._fmax[det].to_value(u.Hz)
            self._freqs[det] = np.linspace(fmin, fmax, self.nstep) * u.Hz
        return self._freqs[det]

    @function_timer
    def bandpass(self, det):
        if det not in self._bandpass:
            # Normalize and interpolate the bandpass
            freqs = self.freqs(det)
            try:
                # If we have a tabulated bandpass, interpolate it
                self._bandpass[det] = np.interp(
                    freqs.to_value(u.Hz),
                    self._bins[det].to_value(u.Hz),
                    self._values[det],
                )
            except AttributeError:
                self._bandpass[det] = np.ones(self.nstep)

            # norm = simpson(self.bandpass[det], x=self.freqs[det])
            norm = integrate_simpson(freqs.to_value(u.Hz), self._bandpass[det])
            if norm == 0:
                raise RuntimeError("Bandpass cannot be normalized")
            self._bandpass[det] /= norm

        return self._bandpass[det]

    @function_timer
    def kcmb2jysr(self, det):
        """Return the unit conversion between K_CMB and Jy/sr"""
        self._get_unit_conversion_coefficients(det)
        return self._kcmb2jysr[det]

    @function_timer
    def kcmb2krj(self, det):
        """Return the unit conversion between K_CMB and K_RJ"""
        self._get_unit_conversion_coefficients(det)
        return self._kcmb2krj[det]

    @function_timer
    def kcmb2w(self, det):
        """Return the unit conversion between K_CMB and W"""
        self._get_unit_conversion_coefficients(det)
        return self._kcmb2w[det]

    @function_timer
    def convolve(self, det, freqs, spectrum, rj=False):
        """Convolve the provided spectrum with the detector bandpass

        Args:
            det(str):  Detector name
            freqs(array of floats):  Spectral bin locations
            spectrum(array of floats):  Spectral bin values
            rj(bool):  Input spectrum is in Rayleigh-Jeans units and
                should be converted into thermal units for convolution

        Returns:
            (array):  The bandpass-convolved spectrum
        """
        freqs_det = self.freqs(det)
        bandpass_det = self.bandpass(det)

        # Interpolate spectrum values to bandpass frequencies
        spectrum_det = np.interp(
            freqs_det.to_value(u.Hz), freqs.to_value(u.Hz), spectrum
        )

        if rj:
            # From brightness to thermodynamic units
            x = h * freqs_det.to_value(u.Hz) / k / TCMB
            rj2cmb = (x / (np.exp(x / 2) - np.exp(-x / 2))) ** -2
            spectrum_det *= rj2cmb

        # Average across the bandpass
        convolved = integrate_simpson(
            freqs_det.to_value(u.Hz), spectrum_det * bandpass_det
        )

        return convolved

_bandpass = {} instance-attribute

_fmax = {} instance-attribute

_fmax_tot = None instance-attribute

_fmin = {} instance-attribute

_fmin_tot = None instance-attribute

_freqs = {} instance-attribute

_kcmb2jysr = {} instance-attribute

_kcmb2krj = {} instance-attribute

_kcmb2w = {} instance-attribute

dets = [] instance-attribute

nstep = nstep instance-attribute

__init__(bandcenters, bandwidths, nstep=101)

All units in GHz

Args

bandcenters(dict) : Dictionary of bandpass centers bandwidths(dict) : Dictionary of bandpass widths nstep(int) : Number of interplation steps to use in convolve()

Source code in toast/instrument.py
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
@function_timer
def __init__(self, bandcenters, bandwidths, nstep=101):
    """All units in GHz

    Args :
        bandcenters(dict) : Dictionary of bandpass centers
        bandwidths(dict) : Dictionary of bandpass widths
        nstep(int) : Number of interplation steps to use in `convolve()`
    """
    self.nstep = nstep
    self.dets = []
    self._fmin = {}
    self._fmax = {}
    for name in bandcenters:
        self.dets.append(name)
        center = bandcenters[name]
        width = bandwidths[name]
        self._fmin[name] = center - 0.5 * width
        self._fmax[name] = center + 0.5 * width
    # The interpolated bandpasses will be cached as needed
    self._fmin_tot = None
    self._fmax_tot = None
    self._freqs = {}
    self._bandpass = {}
    self._kcmb2jysr = {}
    self._kcmb2krj = {}
    self._kcmb2w = {}

_get_unit_conversion_coefficients(det)

Compute and cache the unit conversion coefficients for one detector

Source code in toast/instrument.py
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
@function_timer
def _get_unit_conversion_coefficients(self, det):
    """Compute and cache the unit conversion coefficients for one detector"""

    if (
        det not in self._kcmb2jysr
        or det not in self._kcmb2krj
        or det not in self._kcmb2w
    ):
        # The calculation is a copy from the Hildebrandt and Macias-Perez IDL module for Planck

        nu_cmb = k * TCMB / h
        alpha = 2 * k**3 * TCMB**2 / h**2 / c**2

        cfreq = self.center_frequency(det).to_value(u.Hz)
        freqs = self.freqs(det).to_value(u.Hz)
        bandpass = self.bandpass(det)

        x = freqs / nu_cmb
        db_dt = alpha * x**4 * np.exp(x) / (np.exp(x) - 1) ** 2
        db_dt_rj = 2 * freqs**2 * k / c**2

        self._kcmb2jysr[det] = (
            1e26
            * integrate_simpson(freqs, db_dt * bandpass)
            / integrate_simpson(freqs, cfreq / freqs * bandpass)
        )
        self._kcmb2krj[det] = integrate_simpson(
            freqs, db_dt * bandpass
        ) / integrate_simpson(freqs, db_dt_rj * bandpass)

        # K_CMB->W conversion is from the BoloCalc paper, arXiv:1806.04316
        bandpass = bandpass / np.amax(bandpass)
        self._kcmb2w[det] = integrate_simpson(
            freqs,
            k * (x / (np.exp(x) - 1)) ** 2 * np.exp(x) * bandpass,
        )

    return

bandpass(det)

Source code in toast/instrument.py
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
@function_timer
def bandpass(self, det):
    if det not in self._bandpass:
        # Normalize and interpolate the bandpass
        freqs = self.freqs(det)
        try:
            # If we have a tabulated bandpass, interpolate it
            self._bandpass[det] = np.interp(
                freqs.to_value(u.Hz),
                self._bins[det].to_value(u.Hz),
                self._values[det],
            )
        except AttributeError:
            self._bandpass[det] = np.ones(self.nstep)

        # norm = simpson(self.bandpass[det], x=self.freqs[det])
        norm = integrate_simpson(freqs.to_value(u.Hz), self._bandpass[det])
        if norm == 0:
            raise RuntimeError("Bandpass cannot be normalized")
        self._bandpass[det] /= norm

    return self._bandpass[det]

center_frequency(det, alpha=-1)

Return the effective central frequency for a given spectral index

Source code in toast/instrument.py
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
@function_timer
def center_frequency(self, det, alpha=-1):
    """Return the effective central frequency for a given spectral index"""

    # Which delta function bandpass would produce the same flux density
    freqs = self.freqs(det)
    if alpha == 0:
        # The equation is singular at alpha == 0. Evaluate it on both sides
        # and return the average
        delta = 1e-6
        alpha1 = alpha - delta
        eff1 = self.convolve(det, freqs, freqs.to_value(u.Hz) ** alpha1) ** (
            1 / alpha1
        )
        alpha2 = alpha + delta
        eff2 = self.convolve(det, freqs, freqs.to_value(u.Hz) ** alpha2) ** (
            1 / alpha2
        )
        eff = 0.5 * (eff1 + eff2)
    else:
        # Very simple closed form
        eff = self.convolve(det, freqs, freqs.to_value(u.Hz) ** alpha) ** (
            1 / alpha
        )

    return eff * u.Hz

convolve(det, freqs, spectrum, rj=False)

Convolve the provided spectrum with the detector bandpass

Parameters:

Name Type Description Default
det(str)

Detector name

required
freqs(array of floats

Spectral bin locations

required
spectrum(array of floats

Spectral bin values

required
rj(bool)

Input spectrum is in Rayleigh-Jeans units and should be converted into thermal units for convolution

required

Returns:

Type Description
array

The bandpass-convolved spectrum

Source code in toast/instrument.py
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
@function_timer
def convolve(self, det, freqs, spectrum, rj=False):
    """Convolve the provided spectrum with the detector bandpass

    Args:
        det(str):  Detector name
        freqs(array of floats):  Spectral bin locations
        spectrum(array of floats):  Spectral bin values
        rj(bool):  Input spectrum is in Rayleigh-Jeans units and
            should be converted into thermal units for convolution

    Returns:
        (array):  The bandpass-convolved spectrum
    """
    freqs_det = self.freqs(det)
    bandpass_det = self.bandpass(det)

    # Interpolate spectrum values to bandpass frequencies
    spectrum_det = np.interp(
        freqs_det.to_value(u.Hz), freqs.to_value(u.Hz), spectrum
    )

    if rj:
        # From brightness to thermodynamic units
        x = h * freqs_det.to_value(u.Hz) / k / TCMB
        rj2cmb = (x / (np.exp(x / 2) - np.exp(-x / 2))) ** -2
        spectrum_det *= rj2cmb

    # Average across the bandpass
    convolved = integrate_simpson(
        freqs_det.to_value(u.Hz), spectrum_det * bandpass_det
    )

    return convolved

freqs(det)

Source code in toast/instrument.py
438
439
440
441
442
443
444
@function_timer
def freqs(self, det):
    if det not in self._freqs:
        fmin = self._fmin[det].to_value(u.Hz)
        fmax = self._fmax[det].to_value(u.Hz)
        self._freqs[det] = np.linspace(fmin, fmax, self.nstep) * u.Hz
    return self._freqs[det]

get_range(det=None)

Return the maximum range of frequencies needed for convolution.

Source code in toast/instrument.py
333
334
335
336
337
338
339
340
341
@function_timer
def get_range(self, det=None):
    """Return the maximum range of frequencies needed for convolution."""
    if det is not None:
        return self._fmin[det], self._fmax[det]
    elif self._fmin_tot is None:
        self._fmin_tot = min(self._fmin.values())
        self._fmax_tot = max(self._fmax.values())
    return self._fmin_tot, self._fmax_tot

kcmb2jysr(det)

Return the unit conversion between K_CMB and Jy/sr

Source code in toast/instrument.py
469
470
471
472
473
@function_timer
def kcmb2jysr(self, det):
    """Return the unit conversion between K_CMB and Jy/sr"""
    self._get_unit_conversion_coefficients(det)
    return self._kcmb2jysr[det]

kcmb2krj(det)

Return the unit conversion between K_CMB and K_RJ

Source code in toast/instrument.py
475
476
477
478
479
@function_timer
def kcmb2krj(self, det):
    """Return the unit conversion between K_CMB and K_RJ"""
    self._get_unit_conversion_coefficients(det)
    return self._kcmb2krj[det]

kcmb2w(det)

Return the unit conversion between K_CMB and W

Source code in toast/instrument.py
481
482
483
484
485
@function_timer
def kcmb2w(self, det):
    """Return the unit conversion between K_CMB and W"""
    self._get_unit_conversion_coefficients(det)
    return self._kcmb2w[det]

optical_loading(det, T)

Return the optical loading of a blackbody source. Assumes a diffraction-limited, single-moded polarimeter and perfect optical efficiency arXiv:1806.04316

Parameters:

Name Type Description Default
det(str)

detector name

required
T(float)

Source temperature in Kelvin

required

Returns:

Type Description

(float) : The optical loading in Watts

Source code in toast/instrument.py
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
@function_timer
def optical_loading(self, det, T):
    """Return the optical loading of a blackbody source.
    Assumes a diffraction-limited, single-moded polarimeter
    and perfect optical efficiency
    arXiv:1806.04316

    Args:
        det(str) : detector name
        T(float) : Source temperature in Kelvin

    Returns:
        (float) : The optical loading in Watts
    """

    bandpass = self.bandpass(det)
    # Normalize the bandpass to peak at 1
    bandpass = bandpass / np.amax(bandpass)
    freqs = self.freqs(det).to_value(u.Hz)

    # Power spectral density
    S = h * freqs / (np.exp(h * freqs / k / T) - 1)

    # Integrate over frequency to get power
    power = integrate_simpson(freqs, S * bandpass)

    return power

toast.instrument.Focalplane

Bases: object

Class representing the focalplane for one observation.

The detector_data Table may store arbitrary columns, but several are required. They include:

"name":  The detector name.
"quat":  Each row should be a 4-element numpy array.
"gamma":  If using a half wave plate, we need the rotation angle of the
    detector polarization orientation from the focalplane frame X-axis.

Some columns are optional:

"uid":  Unique integer ID for each detector.  Computed from detector name if
    not specified.
"pol_angle":  Quantity to specify the polarization angle.  Default assumes
    the polarization sensitive direction is aligned with the detector
    quaternion rotation.  Computed if not specified.
"pol_leakage":  Float value "epsilon" between 0-1.  Set to zero by default.
"pol_efficiency":  Float value "eta" = (1 - epsilon) / (1 + epsilon).  Set
    to one by default.
"fwhm":  Quantity with the nominal beam FWHM.  Used for plotting and for
    smoothing of simulated sky signal with PySM.
"bandcenter":  Quantity for the band center.  Used for bandpass integration
    with PySM simulations.
"bandwidth":  Quantity for width of the band.  Used for bandpass integration
    with PySM simulations.
"psd_net":  The detector sensitivity.  Quantity used to create a synthetic
    noise model with the DefaultNoiseModel operator.
"psd_fknee":  Quantity used to create a synthetic noise model with the
    DefaultNoiseModel operator.
"psd_fmin":  Quantity used to create a synthetic noise model with the
    DefaultNoiseModel operator.
"psd_alpha":  Quantity used to create a synthetic noise model with the
    DefaultNoiseModel operator.
"elevation_noise_a" and "elevation_noise_c":  Parameters of elevation scaling
    noise model: PSD_{out} = PSD_{ref} * (a / sin(el) + c)^2.  Only applicable
    to ground data.
"pwv_a0", "pwv_a1" and "pwv_a2":  quadratic fit of the NET modulation by
    PWV.  Only applicable to ground data.

Parameters:

Name Type Description Default
detector_data QTable

Table of detector properties.

None
field_of_view Quantity

Angular diameter of the focal plane. Used to increase the effective size of the focalplane when simulating atmosphere, etc. Will be calculated from the detector offsets by default.

None
sample_rate Quantity

The common (nominal) sample rate for all detectors.

None
thinfp int

Only sample the detectors in the file.

None
Source code in toast/instrument.py
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
class Focalplane(object):
    """Class representing the focalplane for one observation.

    The detector_data Table may store arbitrary columns, but several are required.
    They include:

        "name":  The detector name.
        "quat":  Each row should be a 4-element numpy array.
        "gamma":  If using a half wave plate, we need the rotation angle of the
            detector polarization orientation from the focalplane frame X-axis.

    Some columns are optional:

        "uid":  Unique integer ID for each detector.  Computed from detector name if
            not specified.
        "pol_angle":  Quantity to specify the polarization angle.  Default assumes
            the polarization sensitive direction is aligned with the detector
            quaternion rotation.  Computed if not specified.
        "pol_leakage":  Float value "epsilon" between 0-1.  Set to zero by default.
        "pol_efficiency":  Float value "eta" = (1 - epsilon) / (1 + epsilon).  Set
            to one by default.
        "fwhm":  Quantity with the nominal beam FWHM.  Used for plotting and for
            smoothing of simulated sky signal with PySM.
        "bandcenter":  Quantity for the band center.  Used for bandpass integration
            with PySM simulations.
        "bandwidth":  Quantity for width of the band.  Used for bandpass integration
            with PySM simulations.
        "psd_net":  The detector sensitivity.  Quantity used to create a synthetic
            noise model with the DefaultNoiseModel operator.
        "psd_fknee":  Quantity used to create a synthetic noise model with the
            DefaultNoiseModel operator.
        "psd_fmin":  Quantity used to create a synthetic noise model with the
            DefaultNoiseModel operator.
        "psd_alpha":  Quantity used to create a synthetic noise model with the
            DefaultNoiseModel operator.
        "elevation_noise_a" and "elevation_noise_c":  Parameters of elevation scaling
            noise model: PSD_{out} = PSD_{ref} * (a / sin(el) + c)^2.  Only applicable
            to ground data.
        "pwv_a0", "pwv_a1" and "pwv_a2":  quadratic fit of the NET modulation by
            PWV.  Only applicable to ground data.

    Args:
        detector_data (QTable):  Table of detector properties.
        field_of_view (Quantity):  Angular diameter of the focal plane.  Used to
            increase the effective size of the focalplane when simulating atmosphere,
            etc.  Will be calculated from the detector offsets by default.
        sample_rate (Quantity):  The common (nominal) sample rate for all detectors.
        thinfp (int):  Only sample the detectors in the file.

    """

    XAXIS, YAXIS, ZAXIS = np.eye(3)

    @function_timer
    def __init__(
        self,
        detector_data=None,
        field_of_view=None,
        sample_rate=None,
        thinfp=None,
    ):
        self.detector_data = detector_data
        self.field_of_view = field_of_view
        self.sample_rate = sample_rate
        self.thinfp = thinfp
        if detector_data is not None and len(detector_data) > 0:
            # We have some dets
            self._initialize()

    @function_timer
    def _initialize(self):
        log = Logger.get()

        if self.thinfp is not None:
            # Pick only every `thinfp` pixel on the focal plane
            ndet = len(self.detector_data)
            for idet in range(ndet - 1, -1, -1):
                if int(idet // 2) % self.thinfp != 0:
                    del self.detector_data[idet]

        # Add UID if not given
        if "uid" not in self.detector_data.colnames:
            self.detector_data.add_column(
                Column(
                    name="uid", data=[name_UID(x["name"]) for x in self.detector_data]
                )
            )

        # Build index of detector to table row
        self._det_to_row = {y["name"]: x for x, y in enumerate(self.detector_data)}

        if self.field_of_view is None:
            self._compute_fov()
        self._get_pol_angles()
        self._get_pol_efficiency()
        self._get_bandpass()

    @function_timer
    def _get_bandpass(self):
        """Use the bandpass parameters to instantiate a bandpass model"""

        if "bandcenter" in self.detector_data.colnames:
            bandcenter = {}
            bandwidth = {}
            for row in self.detector_data:
                name = row["name"]
                bandcenter[name] = row["bandcenter"]
                bandwidth[name] = row["bandwidth"]
            self.bandpass = Bandpass(bandcenter, bandwidth)
        else:
            self.bandpass = None
        return

    @function_timer
    def _compute_fov(self):
        """Compute the field of view"""
        # Find the largest distance from the bore sight
        cosangs = list()
        for row in self.detector_data:
            quat = row["quat"]
            vec = qarray.rotate(quat, self.ZAXIS)
            cosangs.append(np.dot(self.ZAXIS, vec))
        mincos = np.amin(cosangs)
        # Add a very small margin to avoid numeric issues
        # in the atmospheric simulation
        self.field_of_view = 1.01 * 2.0 * np.arccos(mincos) * u.radian
        # If we just have boresight detectors, we will need to give this some non-zero
        # value.
        if self.field_of_view == 0:
            self.field_of_view = 1.0 * u.degree

    @function_timer
    def _get_pol_angles(self):
        """Get the detector polarization angles from the quaternions"""

        if "pol_angle" not in self.detector_data.colnames:
            n_rows = len(self.detector_data)
            self.detector_data.add_column(
                Column(name="pol_angle", length=n_rows, unit=u.radian)
            )
            for row in self.detector_data:
                quat = row["quat"]
                a = quat[3]
                d = quat[2]
                pol_angle = np.arctan2(2 * a * d, a**2 - d**2) % np.pi
                row["pol_angle"] = pol_angle * u.radian

    @function_timer
    def _get_pol_efficiency(self):
        """Get the polarization efficiency from polarization leakage or vice versa"""

        n_rows = len(self.detector_data)
        if ("pol_leakage" in self.detector_data.colnames) and (
            "pol_efficiency" in self.detector_data.colnames
        ):
            # Check that efficiency and leakage are consistent
            epsilon = self.detector_data["pol_leakage"]
            eta = self.detector_data["pol_efficiency"]
            np.testing.assert_allclose(
                eta,
                (1 + epsilon) / (1 - epsilon),
                rtol=1e-6,
                err_msg="inconsistent polarization leakage and efficiency",
            )
            return
        elif "pol_leakage" in self.detector_data.colnames:
            self.detector_data.add_column(
                Column(
                    name="pol_efficiency",
                    data=[(1 - x) / (1 + x) for x in self.detector_data["pol_leakage"]],
                )
            )
        elif "pol_efficiency" in self.detector_data.colnames:
            self.detector_data.add_column(
                Column(
                    name="pol_leakage",
                    data=[
                        (1 - x) / (1 + x) for x in self.detector_data["pol_efficiency"]
                    ],
                )
            )
        else:
            self.detector_data.add_column(
                Column(name="pol_efficiency", data=np.ones(n_rows))
            )
            self.detector_data.add_column(
                Column(name="pol_leakage", data=np.zeros(n_rows))
            )

    def __contains__(self, key):
        return key in self._det_to_row

    def __getitem__(self, key):
        return self.detector_data[self._det_to_row[key]]

    def __setitem__(self, key, value):
        if key not in self._det_to_row:
            msg = "cannot assign to non-existent detector '{}'".format(key)
            raise ValueError(msg)
        indx = self._det_to_row[key]
        if hasattr(value, "fields"):
            # numpy structured array
            if value.fields is None:
                raise ValueError("assignment value must be structured")
            for cname, ctype in value.fields.items():
                if cname not in self.detector_data.colnames:
                    msg = "assignment value element '{}' is not a det column".format(
                        cname
                    )
                    raise ValueError(msg)
                self.detector_data[indx][cname] = value[cname]
        elif hasattr(value, "colnames"):
            # table row
            for c in value.colnames:
                if c not in self.detector_data.colnames:
                    msg = "assignment value element '{}' is not a det column".format(c)
                    raise ValueError(msg)
                self.detector_data[indx][c] = value[c]
        else:
            # see if it is like a dictionary
            try:
                for k, v in value.items():
                    if k not in self.detector_data.colnames:
                        msg = (
                            "assignment value element '{}' is not a det column".format(
                                k
                            )
                        )
                        raise ValueError(msg)
                    self.detector_data[indx][k] = v
            except Exception:
                raise ValueError(
                    "assignment value must be a dictionary, Row, or structured array"
                )

    @property
    def detectors(self):
        return list(self._det_to_row.keys())

    @property
    def n_detectors(self):
        return len(self._det_to_row.keys())

    def keys(self):
        return self.detectors

    @function_timer
    def detector_groups(self, column):
        """Group detectors by a common value in one property.

        This returns a dictionary whose keys are the unique values of the specified
        detector_data column.  The values for each key are a list of detectors that
        have that value.  This can be useful for creating detector sets for data
        distribution or for considering detectors with correlations.

        Since the column values will be used for dictionary keys, the column must
        be a data type which is hashable.

        Args:
            column (str):  The detector_data column.

        Returns:
            (dict):  The detector names grouped by unique column values.

        """
        if column not in self.detector_data.colnames:
            raise RuntimeError(f"'{column}' is not a valid det data column")
        detgroups = dict()
        for d in self.detectors:
            indx = self._det_to_row[d]
            val = self.detector_data[column][indx]
            if val not in detgroups:
                detgroups[val] = list()
            detgroups[val].append(d)
        return detgroups

    def __repr__(self):
        value = "<Focalplane: {} detectors, sample_rate = {} Hz, FOV = {} deg, detectors = [".format(
            len(self.detector_data),
            self.sample_rate.to_value(u.Hz),
            self.field_of_view.to_value(u.degree),
        )
        value += "{} .. {}".format(self.detectors[0], self.detectors[-1])
        value += "]>"
        return value

    def __eq__(self, other):
        if self.sample_rate != other.sample_rate:
            return False
        if self.field_of_view != other.field_of_view:
            return False
        if self.detectors != other.detectors:
            return False
        if self.detector_data.colnames != other.detector_data.colnames:
            return False
        if not self.detector_data.values_equal(other.detector_data):
            return False
        return True

    def __ne__(self, other):
        return not self.__eq__(other)

    def load_hdf5(self, handle, comm=None, **kwargs):
        """Load the focalplane from an HDF5 group.

        Args:
            handle (h5py.Group):  The group containing the "focalplane" dataset.
            comm (MPI.Comm):  If loading from a file, optional communicator to broadcast
                across.

        Returns:
            None

        """
        log = Logger.get()

        # Determine if we need to broadcast results.  This occurs if only one process
        # has the file open but the communicator has more than one process.
        need_bcast = hdf5_use_serial(handle, comm)

        if self.detector_data is not None:
            raise RuntimeError("Reading detector data over existing table")

        if handle is not None:
            self.detector_data = read_table_hdf5(handle, path="focalplane")

        if need_bcast and comm is not None:
            self.detector_data = comm.bcast(self.detector_data, root=0)

        # Only use the sampling rate recorded in the file if it was not
        # overridden in the constructor
        if self.sample_rate is None:
            self.sample_rate = self.detector_data.meta["sample_rate"]
        if self.field_of_view is None and "field_of_view" in self.detector_data.meta:
            self.field_of_view = self.detector_data.meta["field_of_view"]

        # Initialize other properties
        self._initialize()

        log.debug_rank(
            f"Focalplane has {len(self.detector_data)} detectors that span "
            f"{self.field_of_view.to_value(u.deg):.3f} degrees and are sampled at "
            f"{self.sample_rate.to_value(u.Hz)} Hz.",
            comm=comm,
        )

    def save_hdf5(self, handle, comm=None, **kwargs):
        """Save the focalplane to an HDF5 group.

        Args:
            handle (h5py.Group):  The parent group of the focalplane dataset.
            comm (MPI.Comm):  If loading from a file, optional communicator to broadcast
                across.

        Returns:
            None

        """
        self.detector_data.meta["sample_rate"] = self.sample_rate
        self.detector_data.meta["field_of_view"] = self.field_of_view
        table_write_parallel_hdf5(self.detector_data, handle, "focalplane", comm=comm)

detector_data = detector_data instance-attribute

detectors property

field_of_view = field_of_view instance-attribute

n_detectors property

sample_rate = sample_rate instance-attribute

thinfp = thinfp instance-attribute

__contains__(key)

Source code in toast/instrument.py
712
713
def __contains__(self, key):
    return key in self._det_to_row

__eq__(other)

Source code in toast/instrument.py
809
810
811
812
813
814
815
816
817
818
819
820
def __eq__(self, other):
    if self.sample_rate != other.sample_rate:
        return False
    if self.field_of_view != other.field_of_view:
        return False
    if self.detectors != other.detectors:
        return False
    if self.detector_data.colnames != other.detector_data.colnames:
        return False
    if not self.detector_data.values_equal(other.detector_data):
        return False
    return True

__getitem__(key)

Source code in toast/instrument.py
715
716
def __getitem__(self, key):
    return self.detector_data[self._det_to_row[key]]

__init__(detector_data=None, field_of_view=None, sample_rate=None, thinfp=None)

Source code in toast/instrument.py
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
@function_timer
def __init__(
    self,
    detector_data=None,
    field_of_view=None,
    sample_rate=None,
    thinfp=None,
):
    self.detector_data = detector_data
    self.field_of_view = field_of_view
    self.sample_rate = sample_rate
    self.thinfp = thinfp
    if detector_data is not None and len(detector_data) > 0:
        # We have some dets
        self._initialize()

__ne__(other)

Source code in toast/instrument.py
822
823
def __ne__(self, other):
    return not self.__eq__(other)

__repr__()

Source code in toast/instrument.py
799
800
801
802
803
804
805
806
807
def __repr__(self):
    value = "<Focalplane: {} detectors, sample_rate = {} Hz, FOV = {} deg, detectors = [".format(
        len(self.detector_data),
        self.sample_rate.to_value(u.Hz),
        self.field_of_view.to_value(u.degree),
    )
    value += "{} .. {}".format(self.detectors[0], self.detectors[-1])
    value += "]>"
    return value

__setitem__(key, value)

Source code in toast/instrument.py
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
def __setitem__(self, key, value):
    if key not in self._det_to_row:
        msg = "cannot assign to non-existent detector '{}'".format(key)
        raise ValueError(msg)
    indx = self._det_to_row[key]
    if hasattr(value, "fields"):
        # numpy structured array
        if value.fields is None:
            raise ValueError("assignment value must be structured")
        for cname, ctype in value.fields.items():
            if cname not in self.detector_data.colnames:
                msg = "assignment value element '{}' is not a det column".format(
                    cname
                )
                raise ValueError(msg)
            self.detector_data[indx][cname] = value[cname]
    elif hasattr(value, "colnames"):
        # table row
        for c in value.colnames:
            if c not in self.detector_data.colnames:
                msg = "assignment value element '{}' is not a det column".format(c)
                raise ValueError(msg)
            self.detector_data[indx][c] = value[c]
    else:
        # see if it is like a dictionary
        try:
            for k, v in value.items():
                if k not in self.detector_data.colnames:
                    msg = (
                        "assignment value element '{}' is not a det column".format(
                            k
                        )
                    )
                    raise ValueError(msg)
                self.detector_data[indx][k] = v
        except Exception:
            raise ValueError(
                "assignment value must be a dictionary, Row, or structured array"
            )

_compute_fov()

Compute the field of view

Source code in toast/instrument.py
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
@function_timer
def _compute_fov(self):
    """Compute the field of view"""
    # Find the largest distance from the bore sight
    cosangs = list()
    for row in self.detector_data:
        quat = row["quat"]
        vec = qarray.rotate(quat, self.ZAXIS)
        cosangs.append(np.dot(self.ZAXIS, vec))
    mincos = np.amin(cosangs)
    # Add a very small margin to avoid numeric issues
    # in the atmospheric simulation
    self.field_of_view = 1.01 * 2.0 * np.arccos(mincos) * u.radian
    # If we just have boresight detectors, we will need to give this some non-zero
    # value.
    if self.field_of_view == 0:
        self.field_of_view = 1.0 * u.degree

_get_bandpass()

Use the bandpass parameters to instantiate a bandpass model

Source code in toast/instrument.py
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
@function_timer
def _get_bandpass(self):
    """Use the bandpass parameters to instantiate a bandpass model"""

    if "bandcenter" in self.detector_data.colnames:
        bandcenter = {}
        bandwidth = {}
        for row in self.detector_data:
            name = row["name"]
            bandcenter[name] = row["bandcenter"]
            bandwidth[name] = row["bandwidth"]
        self.bandpass = Bandpass(bandcenter, bandwidth)
    else:
        self.bandpass = None
    return

_get_pol_angles()

Get the detector polarization angles from the quaternions

Source code in toast/instrument.py
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
@function_timer
def _get_pol_angles(self):
    """Get the detector polarization angles from the quaternions"""

    if "pol_angle" not in self.detector_data.colnames:
        n_rows = len(self.detector_data)
        self.detector_data.add_column(
            Column(name="pol_angle", length=n_rows, unit=u.radian)
        )
        for row in self.detector_data:
            quat = row["quat"]
            a = quat[3]
            d = quat[2]
            pol_angle = np.arctan2(2 * a * d, a**2 - d**2) % np.pi
            row["pol_angle"] = pol_angle * u.radian

_get_pol_efficiency()

Get the polarization efficiency from polarization leakage or vice versa

Source code in toast/instrument.py
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
@function_timer
def _get_pol_efficiency(self):
    """Get the polarization efficiency from polarization leakage or vice versa"""

    n_rows = len(self.detector_data)
    if ("pol_leakage" in self.detector_data.colnames) and (
        "pol_efficiency" in self.detector_data.colnames
    ):
        # Check that efficiency and leakage are consistent
        epsilon = self.detector_data["pol_leakage"]
        eta = self.detector_data["pol_efficiency"]
        np.testing.assert_allclose(
            eta,
            (1 + epsilon) / (1 - epsilon),
            rtol=1e-6,
            err_msg="inconsistent polarization leakage and efficiency",
        )
        return
    elif "pol_leakage" in self.detector_data.colnames:
        self.detector_data.add_column(
            Column(
                name="pol_efficiency",
                data=[(1 - x) / (1 + x) for x in self.detector_data["pol_leakage"]],
            )
        )
    elif "pol_efficiency" in self.detector_data.colnames:
        self.detector_data.add_column(
            Column(
                name="pol_leakage",
                data=[
                    (1 - x) / (1 + x) for x in self.detector_data["pol_efficiency"]
                ],
            )
        )
    else:
        self.detector_data.add_column(
            Column(name="pol_efficiency", data=np.ones(n_rows))
        )
        self.detector_data.add_column(
            Column(name="pol_leakage", data=np.zeros(n_rows))
        )

_initialize()

Source code in toast/instrument.py
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
@function_timer
def _initialize(self):
    log = Logger.get()

    if self.thinfp is not None:
        # Pick only every `thinfp` pixel on the focal plane
        ndet = len(self.detector_data)
        for idet in range(ndet - 1, -1, -1):
            if int(idet // 2) % self.thinfp != 0:
                del self.detector_data[idet]

    # Add UID if not given
    if "uid" not in self.detector_data.colnames:
        self.detector_data.add_column(
            Column(
                name="uid", data=[name_UID(x["name"]) for x in self.detector_data]
            )
        )

    # Build index of detector to table row
    self._det_to_row = {y["name"]: x for x, y in enumerate(self.detector_data)}

    if self.field_of_view is None:
        self._compute_fov()
    self._get_pol_angles()
    self._get_pol_efficiency()
    self._get_bandpass()

detector_groups(column)

Group detectors by a common value in one property.

This returns a dictionary whose keys are the unique values of the specified detector_data column. The values for each key are a list of detectors that have that value. This can be useful for creating detector sets for data distribution or for considering detectors with correlations.

Since the column values will be used for dictionary keys, the column must be a data type which is hashable.

Parameters:

Name Type Description Default
column str

The detector_data column.

required

Returns:

Type Description
dict

The detector names grouped by unique column values.

Source code in toast/instrument.py
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
@function_timer
def detector_groups(self, column):
    """Group detectors by a common value in one property.

    This returns a dictionary whose keys are the unique values of the specified
    detector_data column.  The values for each key are a list of detectors that
    have that value.  This can be useful for creating detector sets for data
    distribution or for considering detectors with correlations.

    Since the column values will be used for dictionary keys, the column must
    be a data type which is hashable.

    Args:
        column (str):  The detector_data column.

    Returns:
        (dict):  The detector names grouped by unique column values.

    """
    if column not in self.detector_data.colnames:
        raise RuntimeError(f"'{column}' is not a valid det data column")
    detgroups = dict()
    for d in self.detectors:
        indx = self._det_to_row[d]
        val = self.detector_data[column][indx]
        if val not in detgroups:
            detgroups[val] = list()
        detgroups[val].append(d)
    return detgroups

keys()

Source code in toast/instrument.py
766
767
def keys(self):
    return self.detectors

load_hdf5(handle, comm=None, **kwargs)

Load the focalplane from an HDF5 group.

Parameters:

Name Type Description Default
handle Group

The group containing the "focalplane" dataset.

required
comm Comm

If loading from a file, optional communicator to broadcast across.

None

Returns:

Type Description

None

Source code in toast/instrument.py
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
def load_hdf5(self, handle, comm=None, **kwargs):
    """Load the focalplane from an HDF5 group.

    Args:
        handle (h5py.Group):  The group containing the "focalplane" dataset.
        comm (MPI.Comm):  If loading from a file, optional communicator to broadcast
            across.

    Returns:
        None

    """
    log = Logger.get()

    # Determine if we need to broadcast results.  This occurs if only one process
    # has the file open but the communicator has more than one process.
    need_bcast = hdf5_use_serial(handle, comm)

    if self.detector_data is not None:
        raise RuntimeError("Reading detector data over existing table")

    if handle is not None:
        self.detector_data = read_table_hdf5(handle, path="focalplane")

    if need_bcast and comm is not None:
        self.detector_data = comm.bcast(self.detector_data, root=0)

    # Only use the sampling rate recorded in the file if it was not
    # overridden in the constructor
    if self.sample_rate is None:
        self.sample_rate = self.detector_data.meta["sample_rate"]
    if self.field_of_view is None and "field_of_view" in self.detector_data.meta:
        self.field_of_view = self.detector_data.meta["field_of_view"]

    # Initialize other properties
    self._initialize()

    log.debug_rank(
        f"Focalplane has {len(self.detector_data)} detectors that span "
        f"{self.field_of_view.to_value(u.deg):.3f} degrees and are sampled at "
        f"{self.sample_rate.to_value(u.Hz)} Hz.",
        comm=comm,
    )

save_hdf5(handle, comm=None, **kwargs)

Save the focalplane to an HDF5 group.

Parameters:

Name Type Description Default
handle Group

The parent group of the focalplane dataset.

required
comm Comm

If loading from a file, optional communicator to broadcast across.

None

Returns:

Type Description

None

Source code in toast/instrument.py
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
def save_hdf5(self, handle, comm=None, **kwargs):
    """Save the focalplane to an HDF5 group.

    Args:
        handle (h5py.Group):  The parent group of the focalplane dataset.
        comm (MPI.Comm):  If loading from a file, optional communicator to broadcast
            across.

    Returns:
        None

    """
    self.detector_data.meta["sample_rate"] = self.sample_rate
    self.detector_data.meta["field_of_view"] = self.field_of_view
    table_write_parallel_hdf5(self.detector_data, handle, "focalplane", comm=comm)

toast.instrument.Telescope

Bases: object

Class representing telescope properties for one observation.

Parameters:

Name Type Description Default
name str

The name of the telescope.

required
uid int

The Unique ID of the telescope. If not specified, constructed from a hash of the site name.

None
focalplane Focalplane

The focalplane for this observation.

None
site Site

The site of the telescope for this observation.

None
Source code in toast/instrument.py
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
class Telescope(object):
    """Class representing telescope properties for one observation.

    Args:
        name (str):  The name of the telescope.
        uid (int):  The Unique ID of the telescope.  If not specified, constructed from
            a hash of the site name.
        focalplane (Focalplane):  The focalplane for this observation.
        site (Site):  The site of the telescope for this observation.

    """

    def __init__(self, name, uid=None, focalplane=None, site=None):
        self.name = name
        self.uid = uid
        if self.uid is None:
            self.uid = name_UID(name)
        if not isinstance(focalplane, Focalplane):
            raise RuntimeError("focalplane should be a Focalplane class instance")
        self.focalplane = focalplane
        if not isinstance(site, Site):
            raise RuntimeError("site should be a Site class instance")
        self.site = site

    def __repr__(self):
        value = "<Telescope '{}': uid = {}, site = {}, ".format(
            self.name,
            self.uid,
            self.site,
        )
        value += "focalplane = {}".format(self.focalplane.__repr__())
        value += ">"
        return value

    def __eq__(self, other):
        if self.name != other.name:
            return False
        if self.uid != other.uid:
            return False
        if self.site != other.site:
            return False
        if self.focalplane != other.focalplane:
            return False
        return True

    def __ne__(self, other):
        return not self.__eq__(other)

focalplane = focalplane instance-attribute

name = name instance-attribute

site = site instance-attribute

uid = uid instance-attribute

__eq__(other)

Source code in toast/instrument.py
979
980
981
982
983
984
985
986
987
988
def __eq__(self, other):
    if self.name != other.name:
        return False
    if self.uid != other.uid:
        return False
    if self.site != other.site:
        return False
    if self.focalplane != other.focalplane:
        return False
    return True

__init__(name, uid=None, focalplane=None, site=None)

Source code in toast/instrument.py
957
958
959
960
961
962
963
964
965
966
967
def __init__(self, name, uid=None, focalplane=None, site=None):
    self.name = name
    self.uid = uid
    if self.uid is None:
        self.uid = name_UID(name)
    if not isinstance(focalplane, Focalplane):
        raise RuntimeError("focalplane should be a Focalplane class instance")
    self.focalplane = focalplane
    if not isinstance(site, Site):
        raise RuntimeError("site should be a Site class instance")
    self.site = site

__ne__(other)

Source code in toast/instrument.py
990
991
def __ne__(self, other):
    return not self.__eq__(other)

__repr__()

Source code in toast/instrument.py
969
970
971
972
973
974
975
976
977
def __repr__(self):
    value = "<Telescope '{}': uid = {}, site = {}, ".format(
        self.name,
        self.uid,
        self.site,
    )
    value += "focalplane = {}".format(self.focalplane.__repr__())
    value += ">"
    return value

toast.instrument.Session

Bases: object

Class representing an observing session.

A session consists of multiple Observation instances with different sets of detectors and possibly different sample rates / times. However these observations are on the same physical telescope and over the same broad time range. A session simply tracks that time range and a unique ID which can be used to group the relevant observations.

Parameters:

Name Type Description Default
name str

The name of the session.

required
uid int

The Unique ID of the session. If not specified, it will be constructed from a hash of the name and the optional start/stop times.

None
start datetime

The overall start of the session.

None
end datetime

The overall end of the session.

None
Source code in toast/instrument.py
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
class Session(object):
    """Class representing an observing session.

    A session consists of multiple Observation instances with different sets of
    detectors and possibly different sample rates / times.  However these
    observations are on the same physical telescope and over the same broad
    time range.  A session simply tracks that time range and a unique ID which
    can be used to group the relevant observations.

    Args:
        name (str):  The name of the session.
        uid (int):  The Unique ID of the session.  If not specified, it will be
            constructed from a hash of the name and the optional start/stop times.
        start (datetime):  The overall start of the session.
        end (datetime):  The overall end of the session.

    """

    def __init__(self, name, uid=None, start=None, end=None):
        self.name = name
        for t in start, end:
            if t is not None and not isinstance(t, datetime.datetime):
                raise RuntimeError("Session start/end must be a datetime or None")
        if uid is not None:
            self.uid = uid
        else:
            # Append start and end times to the session name before
            # evaluating the hash.  This reduces the risk of clashing UIDs.
            session_name = name
            if start is not None:
                session_name += start.ctime()
            if end is not None:
                session_name += end.ctime()
            self.uid = name_UID(session_name)
        self.start = start
        self.end = end

    def __repr__(self):
        value = "<Session '{}': uid = {}, start = {}, end = {}".format(
            self.name, self.uid, self.start, self.end
        )
        value += ">"
        return value

    def __eq__(self, other):
        if self.name != other.name:
            return False
        if self.uid != other.uid:
            return False
        if self.start != other.start:
            return False
        if self.end != other.end:
            return False
        return True

    def __ne__(self, other):
        return not self.__eq__(other)

end = end instance-attribute

name = name instance-attribute

start = start instance-attribute

uid = uid instance-attribute

__eq__(other)

Source code in toast/instrument.py
930
931
932
933
934
935
936
937
938
939
def __eq__(self, other):
    if self.name != other.name:
        return False
    if self.uid != other.uid:
        return False
    if self.start != other.start:
        return False
    if self.end != other.end:
        return False
    return True

__init__(name, uid=None, start=None, end=None)

Source code in toast/instrument.py
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
def __init__(self, name, uid=None, start=None, end=None):
    self.name = name
    for t in start, end:
        if t is not None and not isinstance(t, datetime.datetime):
            raise RuntimeError("Session start/end must be a datetime or None")
    if uid is not None:
        self.uid = uid
    else:
        # Append start and end times to the session name before
        # evaluating the hash.  This reduces the risk of clashing UIDs.
        session_name = name
        if start is not None:
            session_name += start.ctime()
        if end is not None:
            session_name += end.ctime()
        self.uid = name_UID(session_name)
    self.start = start
    self.end = end

__ne__(other)

Source code in toast/instrument.py
941
942
def __ne__(self, other):
    return not self.__eq__(other)

__repr__()

Source code in toast/instrument.py
923
924
925
926
927
928
def __repr__(self):
    value = "<Session '{}': uid = {}, start = {}, end = {}".format(
        self.name, self.uid, self.start, self.end
    )
    value += ">"
    return value