Pointing systematics#
The pointing systematics causes some disturbances in the pointing direction of the detectors. Due to the systematics, we will have a signal from the sky where we are not pointing, and the obtained TOD will be perturbed by the systematics.
These TODs with systematics will be used for map-making, though the pointing matrix in the map-maker will be created with the expected pointings, i.e., pointings without systematics.
In order to simulate the pointing systematics, we multiply the systematic quaternions, which describe perturbations, by the expected quaternions. For example, to add a 5-degree offset around the \(x\)-axis to the positional quaternion \(q_z=(0,0,1,0)\) indicating the \(z\)-axis, the following operation is available:
import litebird_sim as lbs
import numpy as np
def print_quaternion(q):
print("{:.3f} {:.3f} {:.3f} {:.3f}".format(*q))
det = lbs.DetectorInfo(
name="000_000_003_QA_040_T",
sampling_rate_hz=1.0,
)
# Positional quaternion indicating z-axis
q_z = lbs.RotQuaternion(
quats=np.array([0.0, 0.0, 1.0, 0.0]),
)
# Systematic quaternion indicating 5-degree rotation around x-axis
q_sys = lbs.RotQuaternion(
quats=np.array(
lbs.quat_rotation_x(np.deg2rad(5.0))
),
)
# The systematic quaternions are multiplied by the expected quaternions in-place.
lbs.left_multiply_syst_quats(
q_z,
q_sys,
det,
start_time=0.0,
sampling_rate_hz=det.sampling_rate_hz,
)
print("Rotation by 5 deg. around x-axis:")
print_quaternion(q_z.quats[0])
Rotation by 5 deg. around x-axis:
-0.000 -0.044 0.999 -0.000
The function left_multiply_syst_quats() multiplies the systematic quaternions by the expected quaternions in-place.
By using this function, the PointingSys class is implemented to inject the systematic quaternions into the quaternions in a specific coordinate.
PointingSys class#
The PointingSys class is used to inject the systematic quaternions into the quaternions in a specific coordinate.
Note that every coordinate is defined by a left-handed coordinate system, and projected pointings by HEALPix are defined as a view from the spacecraft into the sky.
The supported coordinates are:
FocalplaneCoord: The coordinate of the focal plane with the \(z\)-axis along the boresight.FocalplaneCoord.add_offset(): Add the offset to the detectors in the focal plane.FocalplaneCoord.add_disturb(): Add the time-dependent disturbance to the detectors in the focal plane.
SpacecraftCoord: The coordinate of the spacecraft with the \(z\)-axis along its spin axis.SpacecraftCoord.add_offset(): Add the offset to the entire spacecraft.SpacecraftCoord.add_disturb(): Add the time-dependent disturbance to the entire spacecraft.
HWPCoord: Same coordinates as the focal plane but treats the pointing systematics due to the HWP rotation.HWPCoord.add_hwp_rot_disturb(): Add the rotational disturbance synchronized with the HWP rotation.
These classes are accessible through the PointingSys.focalplane, PointingSys.spacecraft, and PointingSys.hwp attributes.
How to inject the pointing systematics#
It is possible to operate FocalplaneCoord, SpacecraftCoord, and HWPCoord directly, but it is recommended to use the PointingSys class to inject the pointing systematics. The following is an example of injecting vibrations into the focal plane.
import numpy as np
import litebird_sim as lbs
from pathlib import Path
import tomlkit
import healpy as hp
def get_hitmap(nside, pointings):
npix = hp.nside2npix(nside)
ipix = hp.ang2pix(nside, pointings[:, :, 0], pointings[:, :, 1])
hitmap, _ = np.histogram(ipix, bins=np.arange(npix + 1))
return hitmap
start_time = 0.0
sampling_hz = 1.0
telescope = "LFT"
dets = []
# Load the mock focal plane configuration
path_of_toml = (
Path(lbs.__file__).resolve().parent.parent
/ "test"
/ "pointing_sys_reference"
/ "mock_focalplane.toml"
)
with open(path_of_toml, "r", encoding="utf-8") as toml_file:
toml_data = tomlkit.parse(toml_file.read())
for i in range(len(toml_data[telescope])):
dets.append(lbs.DetectorInfo.from_dict(toml_data[telescope][f"det_{i:03}"]))
# Create a simulation object
sim = lbs.Simulation(
base_path="pntsys_example",
start_time=start_time,
duration_s=1000.0,
random_seed=12345,
)
# Set the scanning strategy
sim.set_scanning_strategy(
scanning_strategy=lbs.SpinningScanningStrategy(
spin_sun_angle_rad=np.deg2rad(45.0),
spin_rate_hz=0.05 / 60.0,
precession_rate_hz=1.0 / (3.2 * 60 * 60),
),
delta_time_s=1.0 / sampling_hz,
)
# Set the instrument
sim.set_instrument(
lbs.InstrumentInfo(
name="mock_LiteBIRD",
spin_boresight_angle_rad=np.deg2rad(50.0),
),
)
# Set the HWP
sim.set_hwp(lbs.IdealHWP(sim.instrument.hwp_rpm * 2 * np.pi / 60))
# Create the observations
(obs,) = sim.create_observations(
detectors=dets,
split_list_over_processes=False,
)
# Initialize the pointing systematics object
pntsys = lbs.PointingSys(sim, obs, dets)
nquats = obs.n_samples + 1
axes = ["x", "y"]
noise_sigma_deg = 1.0
for ax in axes:
# Prepare the noise to add vibration to the focal plane
noise_rad_array = np.zeros(nquats)
lbs.add_white_noise(
noise_rad_array, sigma=np.deg2rad(noise_sigma_deg), random=sim.random
)
# Add the vibration to the pointings
pntsys.focalplane.add_disturb(noise_rad_array, ax)
lbs.prepare_pointings(
obs,
sim.instrument,
sim.spin2ecliptic_quats,
sim.hwp,
)
pointings, hwp_angle = obs.get_pointings("all")
nside = 64
hitmap = get_hitmap(nside, pointings)
hp.mollview(hitmap, title="Hit-map with focal plane vibration")
In this code snippet, the pointing systematics are injected into the focal plane by adding white noise with a standard deviation of 1 degree to the focal plane. The hit-map is created with the injected pointing systematics.
Tips for MPI distribution#
In the case we consider time-dependent pointing systematics, we may have to increase the number of quaternions by changing delta_time_s in the Simulation.set_scanning_strategy() method. For example, it is changed by considering the nominal sampling rate of LiteBIRD, i.e., 19 Hz:
sim.set_scanning_strategy(
scanning_strategy=lbs.SpinningScanningStrategy(
spin_sun_angle_rad=np.deg2rad(45.0),
spin_rate_hz=0.05 / 60.0,
precession_rate_hz=1.0 / (3.2 * 60 * 60),
),
delta_time_s=1.0 / 19.0,
)
However, after this execution, all of Simulation.spin2ecliptic_quats from start to end in the observation are calculated and stored in memory. And this is not distributed by MPI, so that even when we use MPI, every process will calculate Simulation.spin2ecliptic_quats from start to end of the observation. This causes a memory issue.
To avoid this issue, a simple trick to distribute the quaternion calculations is:
# After preparing:
# `Simulation` with MPI_COMM_WORLD,
# list of `DetectorInfo`,
# `Imo` object...
# Create observation
(obs,) = sim.create_observations(
detectors=dets,
n_blocks_det=1,
n_blocks_time=size,
split_list_over_processes=False
)
pntsys = lbs.PointingSys(sim, obs, dets)
# Define scanning strategy parameters from IMO
scanning_strategy = lbs.SpinningScanningStrategy.from_imo(
url= f"/releases/vPTEP/satellite/scanning_parameters/",
imo=imo,
)
# Calculate quaternions for each process
sim.spin2ecliptic_quats = scanning_strategy.generate_spin2ecl_quaternions(
start_time=obs.start_time,
time_span_s=obs.n_samples/obs.sampling_rate_hz,
delta_time_s=1.0/19.0,
)
# After injecting the pointing systematics...
lbs.prepare_pointings(
obs,
sim.instrument,
sim.spin2ecliptic_quats,
sim.hwp,
)
Especially, this trick is needed when the time-dependent pointing systematics’ spectrum has higher frequency than the quaternion’s sampling rate given by delta_time_s.
API Reference#
- class litebird_sim.pointing_sys.FocalplaneCoord(sim: Simulation, obs: Observation, detectors: List[DetectorInfo])#
Bases:
objectThis class create an instance of focal plane to add offset and disturbance to the detectors.
Methods in this class multiply systematic quaternions to
Observation.quat(List[RotQuaternion]).- Parameters:
sim (Simulation) –
Simulationinstance.obs (Observation) –
Observationinstance whoseObservation.quatis injected with the systematics.detectors (List[DetectorInfo]) – List of
DetectorInfoinstances.
- add_disturb(noise_rad_matrix: ndarray, axis: str)#
Add a rotational disturbance to the detectors in the focal plane by the specified axis.
This method multiplies systematic quaternions to
Observation.quat(List[RotQuaternion]).If the noise_rad_matrix has the shape [N, t] where N is the number of detectors, t is the number of timestamps, the disturbance will be added to the detectors in the focal plane detector by detector. This represents the independent case where each detector has its own disturbance.
- Parameters:
noise_rad_matrix –
(numpy.ndarray with shape [N, t]): The disturbance to be added to all detectors on the focal plane individually, in the specified direction by axis, in radians.
(numpy.ndarray with shape [, t]): The common-mode disturbance to be added to all detectors on the focal plane in the specified direction by axis, in radians.
axis (str) – The axis in the reference frame around which the rotation is to be performed. It must be one of ‘x’, ‘y’, or ‘z’.
- add_offset(offset_rad, axis: str)#
Add a rotational offset to the detectors in the focal plane by the specified axis.
This method multiplies systematic quaternions to
Observation.quat(List[RotQuaternion]).If the offset_rad is a scalar, it will be added to all the detectors in the focal plane. If the offset_rad is an array with same length of the list of detectors, it will be added to the detectors in the focal plane one by one. In this case, the length of the array must be equal to the number of detectors.
- Parameters:
offset_rad – (in case of float): The same offset to be added to all detectors on the focal plane in the specified direction by axis, in radians. (in case of array): The offset to be added to all dtectors on the focal plane individually, in the specified direction by axis, in radians.
axis (str) – The axis in the reference frame around which the rotation is to be performed. It must be one of ‘x’, ‘y’, or ‘z’.
- class litebird_sim.pointing_sys.HWPCoord(sim: Simulation, obs: Observation, detectors: List[DetectorInfo])#
Bases:
objectA class to add pointing disturbance due to the spinning HWP.
Methods in this class multiply systematic quaternions to
DetectorInfo.quat(RotQuaternion).- Parameters:
sim (Simulation) –
Simulationinstance.obs (Observation) –
Observationinstance whoseObservation.quatis injected with the systematics.detectors (List[DetectorInfo]) – List of
DetectorInfoinstances.
- add_hwp_rot_disturb(tilt_angle_rad: float, ang_speed_radpsec: float, tilt_phase_rad=0.0)#
Add a circular pointing disturbance synchrinized with the HWP to detectors in the focal plane.
After the systematics injection, the pointings will be rotated around an expected pointing direction far from a angular distance of self.tilt_angle_rad. The pointing rotation frequency is determined by self.ang_speed_radpsec.
- Parameters:
tilt_angle_rad (float) – The tilted pointing angle from the expected pointing direction.
ang_speed_radpsec (float) – The angular speed of the pointing circular disturbance.
tilt_phase_rad (float) – Angle at which the gradient direction of HWP’s wedge.
- static get_wedgeHWP_pointing_shift_angle(wedge_angle: float, refractive_index: float)#
Calculate the pointing angle shift due to the wedged HWP. It assumes that the HWP has same refractive index between its ordinary and extraordinary optic axes.
- Parameters:
wedge_angle (float) – angle of the wedge HWP in radian.
refractive_index (float) – refractive index of the HWP.
- Returns:
float, the pointing angle shift due to the wedged HWP.
- class litebird_sim.pointing_sys.PointingSys(sim: Simulation, obs: Observation, detectors: List[DetectorInfo])#
Bases:
objectThis class provide an interface to add offset and disturbance to the instrument and detectors.
- Parameters:
sim (Simulation) –
Simulationinstance whose .spin2ecliptic_quats is injected with the systematics.obs (Observation) –
Observationinstance whoseObservation.quatis injected with the systematics.detectors (List[Detector]) – List of
DetectorInfo.
- class litebird_sim.pointing_sys.SpacecraftCoord(sim: Simulation, obs: Observation, detectors: List[DetectorInfo])#
Bases:
objectThis class create an instance of spacecraft to add offset and disturbance to the instrument.
Methods in this class multiply systematic quaternions to
Simulation.spin2ecliptic_quats(RotQuaternion).- Parameters:
sim (Simulation) –
Simulationinstance.obs (Observation) –
Observationinstance whoseObservation.quatis injected with the systematics.detectors (List[DetectorInfo]) – List of
DetectorInfoinstances.
- add_disturb(noise_rad: ndarray, axis)#
Add a rotational disturbance to the instrument in the spacecraft by the specified axis.
This method multiplies systematic quaternions to
Simulation.spin2ecliptic_quats(RotQuaternion).- Parameters:
noise_rad (1d-numpy.ndarray) – The disturbance to be added in the specified direction by axis, in radians.
axis (str) – The axis in the reference frame around which the rotation is to be performed. It must be one of ‘x’, ‘y’, or ‘z’.
- add_offset(offset_rad, axis: str)#
Add a rotational offset to the instrument in the spacecraft by the specified axis.
This method multiplies systematic quaternions to
Simulation.spin2ecliptic_quats(RotQuaternion).- Parameters:
offset_rad (float) – The offset to be added in the specified direction by axis, in radians.
axis (str) – The axis in the reference frame around which the rotation is to be performed. It must be one of ‘x’, ‘y’, or ‘z’.
- litebird_sim.pointing_sys.get_detector_orientation(detector: DetectorInfo)#
This function returns the orientation of the detector in the focal plane.
- Parameters:
detector (
DetectorInfo) – The detector to get the orientation.
- litebird_sim.pointing_sys.left_multiply_syst_quats(result: RotQuaternion, syst_quats: RotQuaternion, detector: DetectorInfo, start_time, sampling_rate_hz)#
Multiply the systematic quaternions syst_quats to the result quaternions.
- Parameters:
result (
RotQuaternion) – The quaternion to which the rotation is to be added.syst_quats (
RotQuaternion) – The quaternions to be multiplied to the result.detector (
DetectorInfo) – The detector to which the rotation is to be added.start_time – Either a floating-point number or an astropy.time.Time object. It can be None if and only if there is just one quaternion in quats.
sampling_rate_hz – The sampling frequency of the quaternions, in Hertz. It can be None if and only if there is just one quaternion in quats.