Simulated Instrument Signal

In this lesson we cover three TOAST operators for simulating noise and systematics:

  • OpSimNoise for simulating instrumental noise
  • OpSimAtmosphere for simulating atmospheric noise
  • OpSimSSS for simulating scan-synchronous signal (typically ground pick-up) We also introduce TODGround, a synthetic TOD class that simulates constant elevation scanning.
In [1]:
# Are you using a special reservation for a workshop?
# If so, set it here:
nersc_reservation = "toast3"

# Load common tools for all lessons
import sys
sys.path.insert(0, "..")
from lesson_tools import (
    check_nersc,
    fake_focalplane
)
nersc_host, nersc_repo, nersc_resv = check_nersc(reservation=nersc_reservation)

# Capture C++ output in the jupyter cells
%reload_ext wurlitzer
Running on NERSC machine 'cori'
  with access to repos: mp107, desi, m1759, cosmo
Using default repo mp107
Reservation 'toast3' valid from 2019-10-16T09:00:00 to 2019-10-16T17:00:00
Current time is 2019-10-16T05:37:55.852545
Reservation 'toast3' not currently valid

Generating a TOAST data object with one constant elevation scan

We begin by generating an observing schedule containing a single observation. We define one circular patch at (40$^\circ$, -40$^\circ$) RA, DEC with a 10$^\circ$ radius.

In [2]:
! toast_ground_schedule.py \
    --site-lat "-22.958064" \
    --site-lon "-67.786222" \
    --site-alt 5200 \
    --site-name Atacama \
    --telescope LAT \
    --start "2020-01-01 04:00:00" \
    --stop "2020-01-01 05:00:00" \
    --patch-coord C \
    --patch small_patch,1,40,-40,10 \
    --ces-max-time 86400 \
    --out schedule.txt

! cat schedule.txt
TOAST INFO: Adding patch "small_patch"
TOAST INFO: Center-and-width format
TOAST INFO: Global timer: toast_ground_schedule:  0.03 seconds (1 calls)
#Site            Telescope        Latitude [deg] Longitude [deg]   Elevation [m]
 Atacama         LAT                     -22.958         -67.786          5200.0
#Start time UTC       Stop time UTC        Start MJD      Stop MJD       Patch name                          Az min   Az max   El       R/S   Sun el1  Sun az1  Sun el2  Sun az2  Moon el1 Moon az1 Moon el2 Moon az2 Phase Pass  Sub
 2020-01-01 04:00:00  2020-01-01 04:53:00    58849.166667   58849.203472 small_patch                           230.06   242.72    38.03 S       -43.25   190.86   -43.78   174.05    -7.65   256.93   -19.09   251.31  0.32     0   0

For simulating atmospheric noise, we'll also need a TOAST weather file, which contains the historic distributions of temperature, wind, pressure and water vapor at a fixed site. The file is arranged by month and hour of the day to capture seasonal and diurnal variation in these parameters. Here we fetch the weather file for Atacama. There is a another file for Pole.

In [3]:
! [[ ! -e weather_Atacama.fits ]] && wget http://portal.nersc.gov/project/cmb/toast_data/example_data/weather_Atacama.fits
--2019-10-16 05:37:57--  http://portal.nersc.gov/project/cmb/toast_data/example_data/weather_Atacama.fits
Resolving portal.nersc.gov (portal.nersc.gov)... 128.55.201.128
Connecting to portal.nersc.gov (portal.nersc.gov)|128.55.201.128|:80... connected.
HTTP request sent, awaiting response... 301 Moved Permanently
Location: https://portal.nersc.gov/project/cmb/toast_data/example_data/weather_Atacama.fits [following]
--2019-10-16 05:37:57--  https://portal.nersc.gov/project/cmb/toast_data/example_data/weather_Atacama.fits
Connecting to portal.nersc.gov (portal.nersc.gov)|128.55.201.128|:443... connected.
HTTP request sent, awaiting response... 200 OK
Length: 2180160 (2.1M) [image/fits]
Saving to: ‘weather_Atacama.fits’

weather_Atacama.fit 100%[===================>]   2.08M  --.-KB/s    in 0.008s  

2019-10-16 05:37:57 (256 MB/s) - ‘weather_Atacama.fits’ saved [2180160/2180160]

TODGround object

Typical TOAST pipelines like toast_ground_sim.py create synthetic observations by loading the schedule from file and then creating instances of TODGround according to the observing schedule. This is most easily accomplished with toast.pipeline_tools.

In [4]:
import toast
import toast.pipeline_tools
from toast.mpi import MPI

import numpy as np
import matplotlib.pyplot as plt

mpiworld, procs, rank = toast.mpi.get_world()
comm = toast.mpi.Comm(mpiworld)

# A pipeline would create the args object with argparse

class args:
    split_schedule = None
    schedule = "schedule.txt"
    sort_schedule = False  # Matters for parallelization
    weather = "weather_Atacama.fits"
    sample_rate = 10  # Hz
    scan_rate = 1.0  # deg / s
    # Use an artifially low acceleration to show the turn-arounds better
    scan_accel = 0.1  # deg / s^2
    hwp_rpm = None
    hwp_step_deg = None
    hwp_step_time_s = None
    fov = 3.0  # Field-of-view in degrees

# Create a fake focalplane, we could also load one from file.
# The Focalplane class interprets the focalplane dictionary
# created by fake_focalplane() but it can also load the information
# from file.

focalplane = toast.pipeline_tools.Focalplane(
    fake_focalplane(fov=args.fov), sample_rate=args.sample_rate
)

# Load the observing schedule, append weather and focalplane to it
    
schedules = toast.pipeline_tools.load_schedule(args, comm)
toast.pipeline_tools.load_weather(args, comm, schedules)
# There could be more than one observing schedule, but not this time
schedule = schedules[0]
schedule.telescope.focalplane = focalplane

# Create a TODGround object based on the only entry in the schedule

ces = schedule.ceslist[0]  # normally we would loop over entries
totsamples = int((ces.stop_time - ces.start_time) * args.sample_rate)

if comm.comm_group is not None:
    # Available detectors should be split between processes in the group
    ndetrank = comm.comm_group.size
else:
    ndetrank = 1
    
telescope = schedule.telescope  # shorthand
        
tod = toast.todmap.TODGround(
    comm.comm_group,
    telescope.focalplane.detquats,
    totsamples,
    detranks=ndetrank,
    firsttime=ces.start_time,
    rate=args.sample_rate,
    site_lon=telescope.site.lon,
    site_lat=telescope.site.lat,
    site_alt=telescope.site.alt,
    azmin=ces.azmin,
    azmax=ces.azmax,
    el=ces.el,
    scanrate=args.scan_rate,
    scan_accel=args.scan_accel,
    #CES_start=None,
    #CES_stop=None,
    #sun_angle_min=args.sun_angle_min,
    #coord=args.coord,
    #sampsizes=None,
    #report_timing=args.debug,
    hwprpm=args.hwp_rpm,
    hwpstep=args.hwp_step_deg,
    hwpsteptime=args.hwp_step_time_s,
)
TOAST INFO: Load 1 (sub)scans in schedule.txt:  0.00 seconds (1 calls)
TOAST INFO: Loading schedule(s):  0.00 seconds (1 calls)
TOAST INFO: Load weather_Atacama.fits:  0.29 seconds (1 calls)
TOAST INFO: TODGround: simulate scan:  0.01 seconds (1 calls)
TOAST INFO: TODGround: list valid intervals:  0.00 seconds (1 calls)
TOAST INFO: TODGround: call base class constructor:  0.00 seconds (1 calls)
TOAST INFO: TODGround: translate scan pointing:  0.03 seconds (1 calls)

The TODGround objects have capabilities that regular TOD objects do not. Specifically, they can produce detector and boresight pointing in both celestial and horizontal coordinate systems. Here we plot the boresight azimuth and identify the stable science scans.

In [5]:
import matplotlib.pyplot as plt
%matplotlib inline

ind = slice(0, 1000)
times = tod.local_times()[ind]
cflags = tod.local_common_flags()[ind]
az = tod.read_boresight_az()[ind]

turnaround = cflags & tod.TURNAROUND != 0
stable = np.logical_not(turnaround)

plt.plot(times[stable], np.degrees(az[stable]), '.', label="stable")
plt.plot(times[turnaround], np.degrees(az[turnaround]), '.', label="turnaround")
ax = plt.gca()
ax.set_xlabel("UNIX time [s]")
ax.set_ylabel("Azimuth [deg]")
ax.set_title("Boresight azimuth, el = {} deg".format(np.degrees(tod._el)))
plt.legend()
Out[5]:
<matplotlib.legend.Legend at 0x2aaaf411f208>

Here we get the horizontal and celestial pointing and plot the location of every 10th sample

In [6]:
import healpy

healpy.mollview(np.zeros(12), title="Horizontal pointing")
for det in tod.local_dets:
    quat_azel = tod.read_pntg(det, azel=True)[::10]
    az, el = toast.qarray.to_position(quat_azel)
    healpy.projplot(az, el, 'r-')
    healpy.graticule(22.5, verbose=False)

healpy.mollview(np.zeros(12), title="Celestial pointing")
for det in tod.local_dets:
    quat_radec = tod.read_pntg(det)[::10]
    ra, dec = toast.qarray.to_position(quat_radec)
    healpy.projplot(ra, dec, 'r.', ms=1)
    healpy.graticule(22.5, verbose=False)

TOAST data

In [7]:
# Now embed the TOD in an observation dictionary and add other necessary metadata

obs = {}
obs["name"] = "CES-{}-{}-{}-{}-{}".format(
    telescope.site.name, telescope.name, ces.name, ces.scan, ces.subscan
)
obs["tod"] = tod
obs["baselines"] = None
obs["noise"] = telescope.focalplane.noise
obs["id"] = int(ces.mjdstart * 10000)
obs["intervals"] = tod.subscans
obs["site"] = telescope.site
obs["site_name"] = telescope.site.name
obs["site_id"] = telescope.site.id
obs["altitude"] = telescope.site.alt
obs["weather"] = telescope.site.weather
obs["telescope"] = telescope
obs["telescope_name"] = telescope.name
obs["telescope_id"] = telescope.id
obs["focalplane"] = telescope.focalplane.detector_data
obs["fpradius"] = telescope.focalplane.radius
obs["start_time"] = ces.start_time
obs["season"] = ces.season
obs["date"] = ces.start_date
obs["MJD"] = ces.mjdstart
obs["rising"] = ces.rising
obs["mindist_sun"] = ces.mindist_sun
obs["mindist_moon"] = ces.mindist_moon
obs["el_sun"] = ces.el_sun
    
In [8]:
for key, value in obs.items():
    if key == "intervals":
        print("intervals = [{} ... {}] ({} in total)".format(value[0], value[-1], len(value)))
    else:
        print("{} = {}".format(key, value))
name = CES-Atacama-LAT-small_patch-0-0
tod = <TODGround
  14 total detectors and 31800 total samples
  Using MPI communicator None
    In grid dimensions 1 sample ranks x 1 detranks
  Process at (0, 0) in grid has data for:
    Samples 0 - 31799 (inclusive)
    Detectors:
      0A
      0B
      1A
      1B
      2A
      2B
      3A
      3B
      4A
      4B
      5A
      5B
      6A
      6B
    Cache contains 2575800 bytes
>
baselines = None
noise = <toast.tod.sim_noise.AnalyticNoise object at 0x2aaaf51dac50>
id = 588491666
intervals = [<Interval 1577851199.0 - 1577851209.0 (-10 - 89)> ... <Interval 1577854369.1 - 1577854379.2 (31691 - 31791)>] (263 in total)
site = (Site 'Atacama' : ID = 0, lon = -67.786, lat = -22.958, alt = 5200.0 m, weather = (Weather : 'weather_Atacama.fits', site = 0, time = None, year = None, month = None, hour = None, realization = 0))
site_name = Atacama
site_id = 0
altitude = 5200.0
weather = (Weather : 'weather_Atacama.fits', site = 0, time = None, year = None, month = None, hour = None, realization = 0)
telescope = (Telescope 'LAT' : ID = 225, Site = (Site 'Atacama' : ID = 0, lon = -67.786, lat = -22.958, alt = 5200.0 m, weather = (Weather : 'weather_Atacama.fits', site = 0, time = None, year = None, month = None, hour = None, realization = 0)), Focalplane = (Focalplane : 14 detectors, sample_rate = 10 Hz, radius = 1.1558552389176189 deg, detectors = (0A, 0B, 1A, 1B, 2A, 2B, 3A, 3B, 4A, 4B, 5A, 5B, 6A, 6B, ))
telescope_name = LAT
telescope_id = 225
focalplane = {'0A': {'quat': array([0.        , 0.        , 0.38268343, 0.92387953]), 'epsilon': 0, 'rate': 20, 'alpha': 1, 'NET': 1, 'fmin': 0, 'fknee': 0.05, 'fwhm_arcmin': 30}, '0B': {'quat': array([0.        , 0.        , 0.92387953, 0.38268343]), 'epsilon': 0, 'rate': 20, 'alpha': 1, 'NET': 1, 'fmin': 0, 'fknee': 0.05, 'fwhm_arcmin': 30}, '1A': {'quat': array([1.11871541e-18, 1.00764926e-02, 0.00000000e+00, 9.99949231e-01]), 'epsilon': 0, 'rate': 20, 'alpha': 1, 'NET': 1, 'fmin': 0, 'fknee': 0.05, 'fwhm_arcmin': 30}, '1B': {'quat': array([0.00712516, 0.00712516, 0.70707088, 0.70707088]), 'epsilon': 0, 'rate': 20, 'alpha': 1, 'NET': 1, 'fmin': 0, 'fknee': 0.05, 'fwhm_arcmin': 30}, '2A': {'quat': array([-0.00613418,  0.00799422,  0.382664  ,  0.92383263]), 'epsilon': 0, 'rate': 20, 'alpha': 1, 'NET': 1, 'fmin': 0, 'fknee': 0.05, 'fwhm_arcmin': 30}, '2B': {'quat': array([0.00131525, 0.00999029, 0.92383263, 0.382664  ]), 'epsilon': 0, 'rate': 20, 'alpha': 1, 'NET': 1, 'fmin': 0, 'fknee': 0.05, 'fwhm_arcmin': 30}, '3A': {'quat': array([-0.0087265 , -0.00503825,  0.        ,  0.99994923]), 'epsilon': 0, 'rate': 20, 'alpha': 1, 'NET': 1, 'fmin': 0, 'fknee': 0.05, 'fwhm_arcmin': 30}, '3B': {'quat': array([-0.00973314,  0.00260799,  0.70707088,  0.70707088]), 'epsilon': 0, 'rate': 20, 'alpha': 1, 'NET': 1, 'fmin': 0, 'fknee': 0.05, 'fwhm_arcmin': 30}, '4A': {'quat': array([-1.23401444e-18, -1.00764926e-02,  0.00000000e+00,  9.99949231e-01]), 'epsilon': 0, 'rate': 20, 'alpha': 1, 'NET': 1, 'fmin': 0, 'fknee': 0.05, 'fwhm_arcmin': 30}, '4B': {'quat': array([-0.00712516, -0.00712516,  0.70707088,  0.70707088]), 'epsilon': 0, 'rate': 20, 'alpha': 1, 'NET': 1, 'fmin': 0, 'fknee': 0.05, 'fwhm_arcmin': 30}, '5A': {'quat': array([ 0.0087265 , -0.00503825,  0.        ,  0.99994923]), 'epsilon': 0, 'rate': 20, 'alpha': 1, 'NET': 1, 'fmin': 0, 'fknee': 0.05, 'fwhm_arcmin': 30}, '5B': {'quat': array([ 0.00260799, -0.00973314,  0.70707088,  0.70707088]), 'epsilon': 0, 'rate': 20, 'alpha': 1, 'NET': 1, 'fmin': 0, 'fknee': 0.05, 'fwhm_arcmin': 30}, '6A': {'quat': array([0.00999029, 0.00131525, 0.382664  , 0.92383263]), 'epsilon': 0, 'rate': 20, 'alpha': 1, 'NET': 1, 'fmin': 0, 'fknee': 0.05, 'fwhm_arcmin': 30}, '6B': {'quat': array([ 0.00799422, -0.00613418,  0.92383263,  0.382664  ]), 'epsilon': 0, 'rate': 20, 'alpha': 1, 'NET': 1, 'fmin': 0, 'fknee': 0.05, 'fwhm_arcmin': 30}}
fpradius = 1.1558552389176189
start_time = 1577851200.0
season = 2020
date = 2020-01-01
MJD = 58849.166667
rising = False
mindist_sun = 88.71160018273949
mindist_moon = 47.56306383473445
el_sun = -43.25
In [9]:
data = toast.Data(comm)
data.obs.append(obs)

Simulating instrument noise

Instrumental noise simulation is based on the TOAST Noise object appended to each observation.

src/toast/tod/sim_det_noise.py

In [10]:
obs = data.obs[0]
noise = obs["noise"]
print("detectors:", noise.detectors)

# Get and plot noise PSD for det1
det1 = noise.detectors[0]
freq1 = noise.freq(det1)
psd1 = noise.psd(det1)
plt.loglog(freq1, psd1, label=det1)

# get and plot det2 but manipulate the noise PSD to suppress high frequency noise
det2 = noise.detectors[1]
freq2 = noise.freq(det2)
psd2 = noise.psd(det2)
psd2 *= 1e-2 / (freq2 + 1e-2)
plt.loglog(freq2, psd2, label=det2)

ax = plt.gca()
ax.set_xlabel("Frequency [Hz]")
ax.set_ylabel("PSD [K$^2$/Hz]")
ax.set_title("Noise PSD")
plt.legend(loc="best")
detectors: ['0A', '0B', '1A', '1B', '2A', '2B', '3A', '3B', '4A', '4B', '5A', '5B', '6A', '6B']
Out[10]:
<matplotlib.legend.Legend at 0x2aaaf5215f60>
In [11]:
simnoise = toast.tod.OpSimNoise(out="noise", realization=0)
toast.tod.OpCacheClear("noise").exec(data)
simnoise.exec(data)
tod = obs["tod"]
times = tod.local_times()
for det in [det1, det2]:
    simulated = tod.local_signal(det, "noise")
    plt.plot(times, simulated, ".", label=det)
ax.set_xlabel("Unix time [s]")
ax.set_ylabel("Noise [K]")
ax.set_title("Simulated noise")
plt.legend(loc="best")
Out[11]:
<matplotlib.legend.Legend at 0x2aaaf5954a58>

By default the Noise object includes a diagonal mixing matrix and produces uncorrelated detector noise. Let's add a strongly correlated low frequency component to the Noise object.

In [12]:
def print_mixmatrix(noise):
    print("{:10}key".format(""))
    print("{:10}".format("detector"), end="")
    for key in noise.keys:
        print("{:>4}".format(key), end="")
    print()
    print("-" * 80)
    for det in noise.detectors:
        print("{:8} :".format(det), end="")
        for key in noise.keys:
            weight = noise.weight(det, key)
            if np.abs(weight) < 1e-10:
                print("{:4}".format(""), end="")
            else:
                print("{:4}".format(noise.weight(det, key)), end="")
        print()

print("Original mixing matrix")
print_mixmatrix(noise)
Original mixing matrix
          key
detector    0A  0B  1A  1B  2A  2B  3A  3B  4A  4B  5A  5B  6A  6B
--------------------------------------------------------------------------------
0A       :   1                                                    
0B       :       1                                                
1A       :           1                                            
1B       :               1                                        
2A       :                   1                                    
2B       :                       1                                
3A       :                           1                            
3B       :                               1                        
4A       :                                   1                    
4B       :                                       1                
5A       :                                           1            
5B       :                                               1        
6A       :                                                   1    
6B       :                                                       1

We will instantiate a new Noise object with all the original inputs and an additional thermal mode. We could also manipulate the existing noise object but this approach is more representative of how pipelines currently use Noise.

In [13]:
detectors = []
freqs = {}
psds = {}
mixmatrix = {}
indices = {}
for det in noise.detectors:
    detectors.append(det)
    freqs[det] = noise.freq(det).copy()
    psds[det] = noise.psd(det).copy()
    indices[det] = noise.index(det)
    # We do not just copy the old mixing matrix because it is not stored for the trivial case
    mixmatrix[det] = {}
    for key in noise.keys:
        weight = noise.weight(det, key)
        if weight != 0:
            mixmatrix[det][key] = weight

# Create a new PSD and add it to the Noise object inputs
correlated_name = "thermal"
freq = np.logspace(-10, 2)
freq[freq < args.sample_rate / 2]
freq[-1] = args.sample_rate / 2
psd = freq ** -2 * 5e-3
freqs[correlated_name] = freq
psds[correlated_name] = psd
indices[correlated_name] = 999999
for det in noise.detectors:
    mixmatrix[det][correlated_name] = 1

# Create a completely new Noise object
new_noise = toast.tod.Noise(
    detectors=detectors, freqs=freqs, psds=psds, mixmatrix=mixmatrix, indices=indices
)

plt.figure()
for key in [det1, det2, correlated_name]:
    freq = new_noise.freq(key)
    psd = new_noise.psd(key)
    plt.loglog(freq, psd, label=key)
ax = plt.gca()
ax.set_xlabel("Frequency [Hz]")
ax.set_ylabel("PSD [K$^2$/Hz]")
ax.set_title("New noise PSD")
plt.legend(loc="best")
Out[13]:
<matplotlib.legend.Legend at 0x2aaaf5beeac8>
In [14]:
print("\nNew mixing matrix")
print_mixmatrix(new_noise)
New mixing matrix
          key
detector    0A  0B  1A  1B  2A  2B  3A  3B  4A  4B  5A  5B  6A  6Bthermal
--------------------------------------------------------------------------------
0A       :   1                                                       1
0B       :       1                                                   1
1A       :           1                                               1
1B       :               1                                           1
2A       :                   1                                       1
2B       :                       1                                   1
3A       :                           1                               1
3B       :                               1                           1
4A       :                                   1                       1
4B       :                                       1                   1
5A       :                                           1               1
5B       :                                               1           1
6A       :                                                   1       1
6B       :                                                       1   1

We will now replace the Noise in the observation with the one that includes the correlated mode and rerun the simulation

In [15]:
obs["noise"] = new_noise
toast.tod.OpCacheClear("new_noise").exec(data)
simnoise = toast.tod.OpSimNoise(out="new_noise", realization=0)
simnoise.exec(data)

plt.figure()
tod = obs["tod"]
times = tod.local_times()
for det in [det1, det2]:
    simulated = tod.local_signal(det, "new_noise")
    plt.plot(times, simulated, ".", label=det)
ax = plt.gca()
ax.set_xlabel("Unix time [s]")
ax.set_ylabel("Noise [K]")
ax.set_title("Simulated noise with thermal mode")
plt.legend(loc="best")
Out[15]:
<matplotlib.legend.Legend at 0x2aaaf41037f0>

Atmospheric noise simulation

The atmospheric noise module in TOAST

  1. defines a rectangular volume of atmosphere that contains all of the planned lines-of-sight factoring in the wind
  2. compresses the volume by determining which elements are actually needed
  3. builds the element-element covariance matrix matching the Errard. et al. model and instantiates it on the compressed volume
  4. optionally caches the volume elements for future use
  5. simulates the detector signal by performing line-of-sight integrals through the simulated volume while moving the volume with wind

src/toast/todmap/sim_det_atm.py

The atmospheric simulation is based on creating realization of the modified Kolmogorov spectrum (Andrews, L. C. 2004, Field Guide to Atmospheric Optics (SPIE), 112):

In [16]:
dissipation_scale = .01  # in meters
injection_scale = 10.  # in meters
kappamin = 1 / injection_scale
kappamax = 1 / dissipation_scale

n = 1000000
k = np.linspace(1e-2, kappamax * 10, n)

kl = 0.9 / dissipation_scale
k0 = 0.75 / injection_scale

kkl = k / kl
phi = (1 + 1.802 * kkl - .254 * (kkl) ** (7 / 6)) * np.exp(-kkl ** 2) * (k ** 2 + k0 ** 2) ** ( -11 / 6)

plt.figure()
plt.loglog(k, k ** (-11 / 3), "r-", label="$\kappa^{-11/3}$")
plt.loglog(k, phi, "b-", label="Modified")

ax = plt.gca()
for k, label in [(kappamin, "$\kappa_\mathrm{min}$"), (kappamax, "$\kappa_\mathrm{max}$")]:
    ax.axvline(k, linestyle="--", color="black")
    plt.annotate(label, xy=[k * 1.08, 1e-12])

ax.set_xlabel("Wave number [$1/\mathrm{m}$]")
ax.set_ylabel("Power spectrum")
ax.set_title("Modifield Kolmogorov spectrum")

ax.set_xlim([1e-2, 1e3])
ax.set_ylim([1e-13, 1e9])

plt.legend(loc="upper right");

Here we run an example atmospheric simulation for two different injection scales. Be patient, they take about a minute to run.

In [17]:
fig = plt.figure(figsize=[6 * 2, 4 * 3])
name = "atmosphere"

# Delete old atmosphere from disk.  We could also set cachedir=None
! rm -rf atm_cache*

# Run the simulation with two different injection scales and plot the resulting TOD

for iplot, lmax in enumerate([10, 300]):
    cachedir = "atm_cache_{}".format(lmax)
    # Wipe old atmospheric signal
    toast.tod.OpCacheClear(name).exec(data)

    atmsim = toast.todmap.OpSimAtmosphere(
        out=name,
        realization=0,  # Each MC will have a different realization
        zmax=1000,  # maximum altitude to integrate
        lmin_center=0.01,  # Dissipation scale
        lmin_sigma=0,
        lmax_center=lmax,  # Injection scale
        lmax_sigma=0,
        xstep=30, # Volume element size
        ystep=30,
        zstep=30,
        nelem_sim_max=10000,  # Target number of volume elements to consider at a time
        gain=3e-5,  # This gain was calibrated against POLARBEAR data
        # If the wind is strong or the observation is long, the volume becomes
        # too large.  This parameter controls breaking the simulation into
        # disjoint segments
        wind_dist=10000,
        cachedir=cachedir,
        freq=100,
        verbosity=1,  # Print progress to stdout and write out the correlation function
    )
    atmsim.exec(data)
    
    # Pick every second detector because the simulated atmosphere is not polarized
    dets = tod.local_dets[::2]

    # Plot the Kolmogorov correlation function
    ax = fig.add_subplot(3, 2, 1 + iplot)
    kolmo = np.genfromtxt("kolmogorov.txt").T
    ax.semilogx(kolmo[0], kolmo[1])
    ax.set_xlabel("Separation [m]")
    ax.set_ylabel("Correlation factor")
    ax.set_title("Kolmogorov correlation, lmax = {}m".format(lmax))
    ax.grid()

    ax = fig.add_subplot(3, 2, 3 + iplot)
    ind = slice(0, 1000)
    tod = obs["tod"]
    times = tod.local_times()
    for det in dets:
        simulated = tod.local_signal(det, name)
        ax.plot(times[ind], simulated[ind], "-", label=det)
    ax.set_xlabel("Unix time [s]")
    ax.set_ylabel("Atmosphere [K]")
    ax.set_title("Simulated atmosphere, lmax = {}m".format(lmax))
    ax.grid()
    
    ax = fig.add_subplot(3, 2, 5 + iplot)
    for det in dets:
        simulated = tod.local_signal(det, name)
        ax.plot(times, simulated, "-", label=det)
    ax.set_xlabel("Unix time [s]")
    ax.set_ylabel("Atmosphere [K]")
    ax.set_title("Simulated atmosphere, lmax = {}m".format(lmax))
    ax.grid()
    
plt.legend(loc="upper right", bbox_to_anchor=(1.2, 1.00))
plt.subplots_adjust(hspace=0.3)
Creating atm_cache_10/6/6/6
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : Setting up atmosphere simulation
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : Instantiating atmosphere for t = 0.0
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : OpSimAtmosphere: Initialize atmosphere:  0.00 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : Simulating the atmosphere for t = 0.0
TOAST INFO: Kolmogorov initialized in:  1.67 seconds (1 calls)
TOAST INFO: Volume compressed in:  8.44 seconds (1 calls)
27350 / 1330560(2.05553 %) volume elements are needed for the simulation
nx = 140 ny = 88 nz = 108
wx = -1.26998 wy = 0.801063 wz = -0.993287
TOAST INFO: Sparse covariance evaluated in:  0.39 seconds (1 calls)
TOAST INFO: Sparse covariance constructed in:  0.40 seconds (2 calls)
TOAST INFO: Cholesky decomposition done in:  0.82 seconds (1 calls)
 s. N = 10058
TOAST INFO: Realization slice (0 -- 10058) var = 7.52128, constructed in:  0.01 seconds (1 calls)
TOAST INFO: Sparse covariance evaluated in:  0.40 seconds (1 calls)
TOAST INFO: Sparse covariance constructed in:  0.41 seconds (2 calls)
TOAST INFO: Cholesky decomposition done in:  0.31 seconds (1 calls)
 s. N = 10039
TOAST INFO: Realization slice (10058 -- 20097) var = 3.02572, constructed in:  0.01 seconds (1 calls)
TOAST INFO: Sparse covariance evaluated in:  0.20 seconds (1 calls)
TOAST INFO: Sparse covariance constructed in:  0.21 seconds (2 calls)
TOAST INFO: Cholesky decomposition done in:  0.09 seconds (1 calls)
 s. N = 7253
TOAST INFO: Realization slice (20097 -- 27350) var = 1.39633, constructed in:  0.00 seconds (1 calls)
TOAST INFO: Realization constructed in:  2.32 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : OpSimAtmosphere: Simulated atmosphere:  12.44 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : Observing the atmosphere
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : OpSimAtmosphere: Observe atmosphere:  0.15 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : OpSimAtmosphere: Initialize atmosphere:  0.00 seconds (1 calls)
atmsim constructed with 1 processes, 2 threads per process.

Input parameters:
             az = [227.807 - 244.881] (17.0743 degrees)
             el = [36.8741 - 39.1859] (2.31171 degrees)
              t = [1.57785e+09 - 1.57785e+09] (3179.9 s)
           lmin = 0.01 +- 0 m
           lmax = 10 +- 0 m
              w = 1.80032 +- 0 m
           wdir = -172.765 +- 0 degrees 
             z0 = 2000 +- 0 m
             T0 = 276.965 +- 0 K
           zatm = 40000 m
           zmax = 1000 m
       scan frame: 
          xstep = 30 m
          ystep = 30 m
          zstep = 30 m
 horizontal frame: 
         xxstep = 5.14843 m
         yystep = 30 m
         zzstep = 42.1129 m
  nelem_sim_max = 10000
        corrlim = 0.001
      verbosity = 1
           rmin = 0 m
           rmax = 100 m

Atmospheric realization parameters:
 lmin = 0.01 m
 lmax = 10 m
    w = 1.80032 m/s
 easthward wind = -1.78599 m/s
 northward wind = -0.226741 m/s
  az0 = 236.344 degrees
  el0 = 38.03 degrees
   wx = -1.26998 m/s
   wy = 0.801063 m/s
   wz = -0.993287 m/s
 wdir = -172.765 degrees
   z0 = 2000 m
   T0 = 276.965 K

Simulation volume:
   delta_x = 4170.96 m
   delta_y = 2631.66 m
   delta_z = 3223.23 m
Observation cone along the X-axis:
   delta_y_cone = 24.3586 m
   delta_z_cone = 4.68029 m
    xstart = -4038.4 m
    ystart = -42.1793 m
    zstart = -3190.62 m
   maxdist = 102.559 m
        nx = 140
        ny = 88
        nz = 108
        nn = 1330560

Evaluating Kolmogorov correlation at 1000 different separations in range 0 - 5950.59 m
kappamin = 0.1 1/m, kappamax =  100 1/m. nkappa = 1000000
Compressing volume, N = 1330560
Flagged hits, flagging neighbors
Creating compression table
Resizing realization to 27350
X-slice: 0 -- 1560(52 30 m layers) m out of  4200 m indices 0 -- 10058 ( 10058 ) out of 27350

0 : Allocated 3.99903 MB for the sparse covariance matrix. Compression: 0.00518133
0 : Analyzing sparse covariance ... 
0 : Factorizing sparse covariance ... 

0 : Allocated 0.0383682 MB for the sparse factorization.

0 : Allocated 50.8239 MB for the sparse sqrt covariance matrix. Compression: 0.0658498
X-slice: 1560 -- 3000(48 30 m layers) m out of  4200 m indices 10058 -- 20097 ( 10039 ) out of 27350

0 : Allocated 3.98053 MB for the sparse covariance matrix. Compression: 0.00517691
0 : Analyzing sparse covariance ... 
0 : Factorizing sparse covariance ... 

0 : Allocated 0.0382957 MB for the sparse factorization.

0 : Allocated 51.5798 MB for the sparse sqrt covariance matrix. Compression: 0.0670823
X-slice: 3000 -- 4170(39 30 m layers) m out of  4200 m indices 20097 -- 27350 ( 7253 ) out of 27350

0 : Allocated 2.8653 MB for the sparse covariance matrix. Compression: 0.00713912
0 : Analyzing sparse covariance ... 
0 : Factorizing sparse covariance ... 

0 : Allocated 0.027668 MB for the sparse factorization.

0 : Allocated 33.8024 MB for the sparse sqrt covariance matrix. Compression: 0.0842213
Saved metadata to atm_cache_10/6/6/6/14869056_588491666_0_0_metadata.txt
Saved realization to atm_cache_10/6/6/6/14869056_588491666_0_0_realization.dat
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : Simulating the atmosphere for t = 0.0
TOAST INFO: Kolmogorov initialized in:  1.54 seconds (1 calls)
TOAST INFO: Volume compressed in:  0.64 seconds (1 calls)
13171 / 63360(20.7876 %) volume elements are needed for the simulation
nx = 55 ny = 32 nz = 36
wx = -1.26998 wy = 0.801063 wz = -0.993287
TOAST INFO: Sparse covariance evaluated in:  0.37 seconds (1 calls)
TOAST INFO: Sparse covariance constructed in:  0.37 seconds (2 calls)
TOAST INFO: Cholesky decomposition done in:  0.00 seconds (1 calls)
 s. N = 10155
TOAST INFO: Realization slice (0 -- 10155) var = 4.4366, constructed in:  0.00 seconds (1 calls)
TOAST INFO: Sparse covariance evaluated in:  0.03 seconds (1 calls)
TOAST INFO: Sparse covariance constructed in:  0.03 seconds (2 calls)
TOAST INFO: Cholesky decomposition done in:  0.00 seconds (1 calls)
 s. N = 3016
TOAST INFO: Realization slice (10155 -- 13171) var = 1.28794, constructed in:  0.00 seconds (1 calls)
TOAST INFO: Realization constructed in:  0.41 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : OpSimAtmosphere: Simulated atmosphere:  2.60 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : Observing the atmosphere
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
atmsim constructed with 1 processes, 2 threads per process.

Input parameters:
             az = [227.807 - 244.881] (17.0743 degrees)
             el = [36.8741 - 39.1859] (2.31171 degrees)
              t = [1.57785e+09 - 1.57785e+09] (3179.9 s)
           lmin = 0.01 +- 0 m
           lmax = 10 +- 0 m
              w = 1.80032 +- 0 m
           wdir = -172.765 +- 0 degrees 
             z0 = 2000 +- 0 m
             T0 = 276.965 +- 0 K
           zatm = 40000 m
           zmax = 1000 m
       scan frame: 
          xstep = 94.8683 m
          ystep = 94.8683 m
          zstep = 94.8683 m
 horizontal frame: 
         xxstep = 16.2808 m
         yystep = 94.8683 m
         zzstep = 133.173 m
  nelem_sim_max = 10000
        corrlim = 0.001
      verbosity = 1
           rmin = 100 m
           rmax = 1000 m

Atmospheric realization parameters:
 lmin = 0.01 m
 lmax = 10 m
    w = 1.80032 m/s
 easthward wind = -1.78599 m/s
 northward wind = -0.226741 m/s
  az0 = 236.344 degrees
  el0 = 38.03 degrees
   wx = -1.26998 m/s
   wy = 0.801063 m/s
   wz = -0.993287 m/s
 wdir = -172.765 degrees
   z0 = 2000 m
   T0 = 276.965 K

Simulation volume:
   delta_x = 5158.86 m
   delta_y = 2980.62 m
   delta_z = 3395.09 m
Observation cone along the X-axis:
   delta_y_cone = 243.586 m
   delta_z_cone = 46.8029 m
    xstart = -4038.4 m
    ystart = -216.661 m
    zstart = -3274.11 m
   maxdist = 1025.59 m
        nx = 55
        ny = 32
        nz = 36
        nn = 63360

Evaluating Kolmogorov correlation at 1000 different separations in range 0 - 6926.02 m
kappamin = 0.1 1/m, kappamax =  100 1/m. nkappa = 1000000
Compressing volume, N = 63360
Flagged hits, flagging neighbors
Creating compression table
Resizing realization to 13171
X-slice: 0 -- 3605(38 94.8683 m layers) m out of  5217.76 m indices 0 -- 10155 ( 10155 ) out of 13171

0 : Allocated 0.154953 MB for the sparse covariance matrix. Compression: 0.000196947
0 : Analyzing sparse covariance ... 
0 : Factorizing sparse covariance ... 

0 : Allocated 0.154953 MB for the sparse factorization.

0 : Allocated 0.154953 MB for the sparse sqrt covariance matrix. Compression: 0.000196947
X-slice: 3605 -- 5122.89(16 94.8683 m layers) m out of  5217.76 m indices 10155 -- 13171 ( 3016 ) out of 13171

0 : Allocated 0.0460205 MB for the sparse covariance matrix. Compression: 0.00066313
0 : Analyzing sparse covariance ... 
0 : Factorizing sparse covariance ... 

0 : Allocated 0.0460205 MB for the sparse factorization.

0 : Allocated 0.0460205 MB for the sparse sqrt covariance matrix. Compression: 0.00066313
Saved metadata to atm_cache_10/6/6/6/14869056_588491666_1_0_metadata.txt
Saved realization to atm_cache_10/6/6/6/14869056_588491666_1_0_realization.dat
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : OpSimAtmosphere: Observe atmosphere:  0.15 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : OpSimAtmosphere: Initialize atmosphere:  0.00 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : Simulating the atmosphere for t = 0.0
TOAST INFO: Kolmogorov initialized in:  1.66 seconds (1 calls)
TOAST INFO: Volume compressed in:  0.03 seconds (1 calls)
2364 / 3120(75.7692 %) volume elements are needed for the simulation
nx = 20 ny = 12 nz = 13
wx = -1.26998 wy = 0.801063 wz = -0.993287
TOAST INFO: Sparse covariance evaluated in:  0.02 seconds (1 calls)
TOAST INFO: Sparse covariance constructed in:  0.02 seconds (2 calls)
TOAST INFO: Cholesky decomposition done in:  0.00 seconds (1 calls)
 s. N = 2364
TOAST INFO: Realization slice (0 -- 2364) var = 3.42455, constructed in:  0.00 seconds (1 calls)
TOAST INFO: Realization constructed in:  0.02 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : OpSimAtmosphere: Simulated atmosphere:  1.72 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : Observing the atmosphere
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : OpSimAtmosphere: Observe atmosphere:  0.11 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : Simulated and observed atmosphere:  17.16 seconds (1 calls)
Creating atm_cache_300/6/6/6
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : Setting up atmosphere simulation
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : Instantiating atmosphere for t = 0.0
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : OpSimAtmosphere: Initialize atmosphere:  0.00 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : Simulating the atmosphere for t = 0.0
TOAST INFO: Kolmogorov initialized in:  1.55 seconds (1 calls)
TOAST INFO: Volume compressed in:  7.96 seconds (1 calls)
27350 / 1330560(2.05553 %) volume elements are needed for the simulation
nx = 140 ny = 88 nz = 108
wx = -1.26998 wy = 0.801063 wz = -0.993287
TOAST INFO: Sparse covariance evaluated in:  7.99 seconds (1 calls)
TOAST INFO: Sparse covariance constructed in:  12.39 seconds (2 calls)
TOAST INFO: Cholesky decomposition done in:  14.49 seconds (1 calls)
 s. N = 10058
TOAST INFO: Realization slice (0 -- 10058) var = 4.44867, constructed in:  0.07 seconds (1 calls)
TOAST INFO: Sparse covariance evaluated in:  4.06 seconds (1 calls)
TOAST INFO: Sparse covariance constructed in:  8.32 seconds (2 calls)
TOAST INFO: Cholesky decomposition done in:  15.67 seconds (1 calls)
 s. N = 10039
TOAST INFO: Realization slice (10058 -- 20097) var = 2.08106, constructed in:  0.09 seconds (1 calls)
TOAST INFO: Sparse covariance evaluated in:  2.12 seconds (1 calls)
TOAST INFO: Sparse covariance constructed in:  4.54 seconds (2 calls)
TOAST INFO: Cholesky decomposition done in:  7.07 seconds (1 calls)
 s. N = 7253
TOAST INFO: Realization slice (20097 -- 27350) var = 0.892358, constructed in:  0.04 seconds (1 calls)
TOAST INFO: Realization constructed in:  63.34 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : OpSimAtmosphere: Simulated atmosphere:  72.86 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : Observing the atmosphere
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : OpSimAtmosphere: Observe atmosphere:  0.13 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : OpSimAtmosphere: Initialize atmosphere:  0.00 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : Simulating the atmosphere for t = 0.0
TOAST INFO: Kolmogorov initialized in:  1.55 seconds (1 calls)
TOAST INFO: Volume compressed in:  0.59 seconds (1 calls)
13171 / 63360(20.7876 %) volume elements are needed for the simulation
nx = 55 ny = 32 nz = 36
wx = -1.26998 wy = 0.801063 wz = -0.993287
TOAST INFO: Sparse covariance evaluated in:  6.84 seconds (1 calls)
TOAST INFO: Sparse covariance constructed in:  10.61 seconds (2 calls)
TOAST INFO: Cholesky decomposition done in:  13.14 seconds (1 calls)
 s. N = 10155
TOAST INFO: Realization slice (0 -- 10155) var = 4.11504, constructed in:  0.09 seconds (1 calls)
TOAST INFO: Sparse covariance evaluated in:  0.32 seconds (1 calls)
TOAST INFO: Sparse covariance constructed in:  0.64 seconds (2 calls)
TOAST INFO: Cholesky decomposition done in:  0.78 seconds (1 calls)
 s. N = 3016
TOAST INFO: Realization slice (10155 -- 13171) var = 1.10595, constructed in:  0.01 seconds (1 calls)
TOAST INFO: Realization constructed 
atmsim constructed with 1 processes, 2 threads per process.

Input parameters:
             az = [227.807 - 244.881] (17.0743 degrees)
             el = [36.8741 - 39.1859] (2.31171 degrees)
              t = [1.57785e+09 - 1.57785e+09] (3179.9 s)
           lmin = 0.01 +- 0 m
           lmax = 10 +- 0 m
              w = 1.80032 +- 0 m
           wdir = -172.765 +- 0 degrees 
             z0 = 2000 +- 0 m
             T0 = 276.965 +- 0 K
           zatm = 40000 m
           zmax = 1000 m
       scan frame: 
          xstep = 300 m
          ystep = 300 m
          zstep = 300 m
 horizontal frame: 
         xxstep = 51.4843 m
         yystep = 300 m
         zzstep = 421.129 m
  nelem_sim_max = 10000
        corrlim = 0.001
      verbosity = 1
           rmin = 1000 m
           rmax = 10000 m

Atmospheric realization parameters:
 lmin = 0.01 m
 lmax = 10 m
    w = 1.80032 m/s
 easthward wind = -1.78599 m/s
 northward wind = -0.226741 m/s
  az0 = 236.344 degrees
  el0 = 38.03 degrees
   wx = -1.26998 m/s
   wy = 0.801063 m/s
   wz = -0.993287 m/s
 wdir = -172.765 degrees
   z0 = 2000 m
   T0 = 276.965 K

Simulation volume:
   delta_x = 5961.58 m
   delta_y = 3532.82 m
   delta_z = 3832.63 m
Observation cone along the X-axis:
   delta_y_cone = 385.519 m
   delta_z_cone = 74.0743 m
    xstart = -4038.4 m
    ystart = -492.76 m
    zstart = -3491.29 m
   maxdist = 1623.18 m
        nx = 20
        ny = 12
        nz = 13
        nn = 3120

Evaluating Kolmogorov correlation at 1000 different separations in range 0 - 7998.17 m
kappamin = 0.1 1/m, kappamax =  100 1/m. nkappa = 1000000
Compressing volume, N = 3120
Flagged hits, flagging neighbors
Creating compression table
Resizing realization to 2364
X-slice: 0 -- 5700(19 300 m layers) m out of  6000 m indices 0 -- 2364 ( 2364 ) out of 2364

0 : Allocated 0.0360718 MB for the sparse covariance matrix. Compression: 0.000846024
0 : Analyzing sparse covariance ... 
0 : Factorizing sparse covariance ... 

0 : Allocated 0.0360718 MB for the sparse factorization.

0 : Allocated 0.0360718 MB for the sparse sqrt covariance matrix. Compression: 0.000846024
Saved metadata to atm_cache_10/6/6/6/14869056_588491666_2_0_metadata.txt
Saved realization to atm_cache_10/6/6/6/14869056_588491666_2_0_realization.dat
atmsim constructed with 1 processes, 2 threads per process.

Input parameters:
             az = [227.807 - 244.881] (17.0743 degrees)
             el = [36.8741 - 39.1859] (2.31171 degrees)
              t = [1.57785e+09 - 1.57785e+09] (3179.9 s)
           lmin = 0.01 +- 0 m
           lmax = 300 +- 0 m
              w = 1.80032 +- 0 m
           wdir = -172.765 +- 0 degrees 
             z0 = 2000 +- 0 m
             T0 = 276.965 +- 0 K
           zatm = 40000 m
           zmax = 1000 m
       scan frame: 
          xstep = 30 m
          ystep = 30 m
          zstep = 30 m
 horizontal frame: 
         xxstep = 5.14843 m
         yystep = 30 m
         zzstep = 42.1129 m
  nelem_sim_max = 10000
        corrlim = 0.001
      verbosity = 1
           rmin = 0 m
           rmax = 100 m

Atmospheric realization parameters:
 lmin = 0.01 m
 lmax = 300 m
    w = 1.80032 m/s
 easthward wind = -1.78599 m/s
 northward wind = -0.226741 m/s
  az0 = 236.344 degrees
  el0 = 38.03 degrees
   wx = -1.26998 m/s
   wy = 0.801063 m/s
   wz = -0.993287 m/s
 wdir = -172.765 degrees
   z0 = 2000 m
   T0 = 276.965 K

Simulation volume:
   delta_x = 4170.96 m
   delta_y = 2631.66 m
   delta_z = 3223.23 m
Observation cone along the X-axis:
   delta_y_cone = 24.3586 m
   delta_z_cone = 4.68029 m
    xstart = -4038.4 m
    ystart = -42.1793 m
    zstart = -3190.62 m
   maxdist = 102.559 m
        nx = 140
        ny = 88
        nz = 108
        nn = 1330560

Evaluating Kolmogorov correlation at 1000 different separations in range 0 - 5950.59 m
kappamin = 0.00333333 1/m, kappamax =  100 1/m. nkappa = 1000000
Compressing volume, N = 1330560
Flagged hits, flagging neighbors
Creating compression table
Resizing realization to 27350
X-slice: 0 -- 1560(52 30 m layers) m out of  4200 m indices 0 -- 10058 ( 10058 ) out of 27350

0 : Allocated 578.957 MB for the sparse covariance matrix. Compression: 0.750124
0 : Analyzing sparse covariance ... 
0 : Factorizing sparse covariance ... 

0 : Allocated 0.0383682 MB for the sparse factorization.

0 : Allocated 578.957 MB for the sparse sqrt covariance matrix. Compression: 0.750124
X-slice: 1560 -- 3000(48 30 m layers) m out of  4200 m indices 10058 -- 20097 ( 10039 ) out of 27350

0 : Allocated 576.772 MB for the sparse covariance matrix. Compression: 0.750125
0 : Analyzing sparse covariance ... 
0 : Factorizing sparse covariance ... 

0 : Allocated 0.0382957 MB for the sparse factorization.

0 : Allocated 576.772 MB for the sparse sqrt covariance matrix. Compression: 0.750125
X-slice: 3000 -- 4170(39 30 m layers) m out of  4200 m indices 20097 -- 27350 ( 7253 ) out of 27350

0 : Allocated 301.083 MB for the sparse covariance matrix. Compression: 0.750172
0 : Analyzing sparse covariance ... 
0 : Factorizing sparse covariance ... 

0 : Allocated 0.027668 MB for the sparse factorization.

0 : Allocated 301.083 MB for the sparse sqrt covariance matrix. Compression: 0.750172
Saved metadata to atm_cache_300/6/6/6/14869056_588491666_0_0_metadata.txt
Saved realization to atm_cache_300/6/6/6/14869056_588491666_0_0_realization.dat
in:  25.59 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : OpSimAtmosphere: Simulated atmosphere:  27.74 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : Observing the atmosphere
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
atmsim constructed with 1 processes, 2 threads per process.

Input parameters:
             az = [227.807 - 244.881] (17.0743 degrees)
             el = [36.8741 - 39.1859] (2.31171 degrees)
              t = [1.57785e+09 - 1.57785e+09] (3179.9 s)
           lmin = 0.01 +- 0 m
           lmax = 300 +- 0 m
              w = 1.80032 +- 0 m
           wdir = -172.765 +- 0 degrees 
             z0 = 2000 +- 0 m
             T0 = 276.965 +- 0 K
           zatm = 40000 m
           zmax = 1000 m
       scan frame: 
          xstep = 94.8683 m
          ystep = 94.8683 m
          zstep = 94.8683 m
 horizontal frame: 
         xxstep = 16.2808 m
         yystep = 94.8683 m
         zzstep = 133.173 m
  nelem_sim_max = 10000
        corrlim = 0.001
      verbosity = 1
           rmin = 100 m
           rmax = 1000 m

Atmospheric realization parameters:
 lmin = 0.01 m
 lmax = 300 m
    w = 1.80032 m/s
 easthward wind = -1.78599 m/s
 northward wind = -0.226741 m/s
  az0 = 236.344 degrees
  el0 = 38.03 degrees
   wx = -1.26998 m/s
   wy = 0.801063 m/s
   wz = -0.993287 m/s
 wdir = -172.765 degrees
   z0 = 2000 m
   T0 = 276.965 K

Simulation volume:
   delta_x = 5158.86 m
   delta_y = 2980.62 m
   delta_z = 3395.09 m
Observation cone along the X-axis:
   delta_y_cone = 243.586 m
   delta_z_cone = 46.8029 m
    xstart = -4038.4 m
    ystart = -216.661 m
    zstart = -3274.11 m
   maxdist = 1025.59 m
        nx = 55
        ny = 32
        nz = 36
        nn = 63360

Evaluating Kolmogorov correlation at 1000 different separations in range 0 - 6926.02 m
kappamin = 0.00333333 1/m, kappamax =  100 1/m. nkappa = 1000000
Compressing volume, N = 63360
Flagged hits, flagging neighbors
Creating compression table
Resizing realization to 13171
X-slice: 0 -- 3605(38 94.8683 m layers) m out of  5217.76 m indices 0 -- 10155 ( 10155 ) out of 13171

0 : Allocated 477.712 MB for the sparse covariance matrix. Compression: 0.607179
0 : Analyzing sparse covariance ... 
0 : Factorizing sparse covariance ... 

0 : Allocated 0.0387383 MB for the sparse factorization.

0 : Allocated 590.177 MB for the sparse sqrt covariance matrix. Compression: 0.750123
X-slice: 3605 -- 5122.89(16 94.8683 m layers) m out of  5217.76 m indices 10155 -- 13171 ( 3016 ) out of 13171

0 : Allocated 52.0632 MB for the sparse covariance matrix. Compression: 0.750202
0 : Analyzing sparse covariance ... 
0 : Factorizing sparse covariance ... 

0 : Allocated 0.0115051 MB for the sparse factorization.

0 : Allocated 52.078 MB for the sparse sqrt covariance matrix. Compression: 0.750414
Saved metadata to atm_cache_300/6/6/6/14869056_588491666_1_0_metadata.txt
Saved realization to atm_cache_300/6/6/6/14869056_588491666_1_0_realization.dat
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO:  samples observed in:  0.01 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : OpSimAtmosphere: Observe atmosphere:  0.24 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : OpSimAtmosphere: Initialize atmosphere:  0.00 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : Simulating the atmosphere for t = 0.0
TOAST INFO: Kolmogorov initialized in:  1.49 seconds (1 calls)
TOAST INFO: Volume compressed in:  0.03 seconds (1 calls)
2364 / 3120(75.7692 %) volume elements are needed for the simulation
nx = 20 ny = 12 nz = 13
wx = -1.26998 wy = 0.801063 wz = -0.993287
TOAST INFO: Sparse covariance evaluated in:  0.30 seconds (1 calls)
TOAST INFO: Sparse covariance constructed in:  0.38 seconds (2 calls)
TOAST INFO: Cholesky decomposition done in:  0.23 seconds (1 calls)
 s. N = 2364
TOAST INFO: Realization slice (0 -- 2364) var = 3.07798, constructed in:  0.00 seconds (1 calls)
TOAST INFO: Realization constructed in:  0.63 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : OpSimAtmosphere: Simulated atmosphere:  2.16 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : Observing the atmosphere
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO:  samples observed in:  0.00 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : OpSimAtmosphere: Observe atmosphere:  0.19 seconds (1 calls)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : Simulated and observed atmosphere:  103.33 seconds (1 calls)
atmsim constructed with 1 processes, 2 threads per process.

Input parameters:
             az = [227.807 - 244.881] (17.0743 degrees)
             el = [36.8741 - 39.1859] (2.31171 degrees)
              t = [1.57785e+09 - 1.57785e+09] (3179.9 s)
           lmin = 0.01 +- 0 m
           lmax = 300 +- 0 m
              w = 1.80032 +- 0 m
           wdir = -172.765 +- 0 degrees 
             z0 = 2000 +- 0 m
             T0 = 276.965 +- 0 K
           zatm = 40000 m
           zmax = 1000 m
       scan frame: 
          xstep = 300 m
          ystep = 300 m
          zstep = 300 m
 horizontal frame: 
         xxstep = 51.4843 m
         yystep = 300 m
         zzstep = 421.129 m
  nelem_sim_max = 10000
        corrlim = 0.001
      verbosity = 1
           rmin = 1000 m
           rmax = 10000 m

Atmospheric realization parameters:
 lmin = 0.01 m
 lmax = 300 m
    w = 1.80032 m/s
 easthward wind = -1.78599 m/s
 northward wind = -0.226741 m/s
  az0 = 236.344 degrees
  el0 = 38.03 degrees
   wx = -1.26998 m/s
   wy = 0.801063 m/s
   wz = -0.993287 m/s
 wdir = -172.765 degrees
   z0 = 2000 m
   T0 = 276.965 K

Simulation volume:
   delta_x = 5961.58 m
   delta_y = 3532.82 m
   delta_z = 3832.63 m
Observation cone along the X-axis:
   delta_y_cone = 385.519 m
   delta_z_cone = 74.0743 m
    xstart = -4038.4 m
    ystart = -492.76 m
    zstart = -3491.29 m
   maxdist = 1623.18 m
        nx = 20
        ny = 12
        nz = 13
        nn = 3120

Evaluating Kolmogorov correlation at 1000 different separations in range 0 - 7998.17 m
kappamin = 0.00333333 1/m, kappamax =  100 1/m. nkappa = 1000000
Compressing volume, N = 3120
Flagged hits, flagging neighbors
Creating compression table
Resizing realization to 2364
X-slice: 0 -- 5700(19 300 m layers) m out of  6000 m indices 0 -- 2364 ( 2364 ) out of 2364

0 : Allocated 18.829 MB for the sparse covariance matrix. Compression: 0.441613
0 : Analyzing sparse covariance ... 
0 : Factorizing sparse covariance ... 

0 : Allocated 0.00901794 MB for the sparse factorization.

0 : Allocated 32.0002 MB for the sparse sqrt covariance matrix. Compression: 0.750529
Saved metadata to atm_cache_300/6/6/6/14869056_588491666_2_0_metadata.txt
Saved realization to atm_cache_300/6/6/6/14869056_588491666_2_0_realization.dat

Ground (scan-synchronous) signal

TOAST includes a new module for simulating scan-synchronous signals. It works by sampling a provided or on-the-fly synthesized low resolution map using the horizontal detector pointing instead of sky pointing.

src/toast/todmap/sss.py

In [18]:
toast.tod.OpCacheClear("sss").exec(data)

sim_sss = toast.todmap.OpSimScanSynchronousSignal(
    out="sss",
    realization=0,
    nside=512,    # internal resolution for the simulated ground map
    fwhm=3,      # smoothing scale for the ground map
    scale=1e-3,   # RMS for observations at el=45 deg
    lmax=512,     # expansion lmax
    power=-1,     # power law that suppresses ground signal at higher elevation
    path=None,    # alternative ground map in Healpix format, will override nside, scale, lmax and power
)

sim_sss.exec(data)

tod = obs["tod"]
times = tod.local_times()
ind = slice(0, 600)
for det in tod.local_dets[::2]:
    simulated = tod.local_signal(det, "sss")
    plt.plot(times[ind], simulated[ind], ".", label=det)
ax = plt.gca()
ax.set_xlabel("Unix time [s]")
ax.set_ylabel("Atmosphere [K]")
ax.set_title("Simulated ground signal")
plt.legend(loc="best")

# The ground map is simulated on-the-fly and discarded afterwards.  By default, each observation gets
# an entirely independent ground map.  Here is an example:

weather = obs["weather"]
weather.set(0, 0, 0)
key1, key2, counter1, counter2 = sim_sss._get_rng_keys(obs)
ground_map = sim_sss._simulate_sss(key1, key2, counter1, counter2, weather, mpiworld)

healpy.mollview(ground_map, cmap="coolwarm", title="Simulated ground map + hits")
for det in tod.local_dets[::2]:
    quat_azel = tod.read_pntg(det, azel=True)
    theta, phi = toast.qarray.to_position(quat_azel)
    healpy.projplot(theta, phi, '-')
healpy.graticule(22.5, verbose=False)

healpy.gnomview(
    ground_map,
    cmap="coolwarm",
    title="Simulated ground map + hits",
    rot=[360 - .5 * (ces.azmin + ces.azmax), ces.el],
    reso=8,
)
for det in tod.local_dets[::2]:
    quat_azel = tod.read_pntg(det, azel=True)
    theta, phi = toast.qarray.to_position(quat_azel)
    healpy.projplot(theta, phi, 'k.')
healpy.graticule(5, verbose=False)
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : Setting up SSS simulation
Sigma is 76.438962 arcmin (0.022235 rad) 
-> fwhm is 180.000000 arcmin
Sigma is 0.000000 arcmin (0.000000 rad) 
-> fwhm is 0.000000 arcmin
TOAST INFO: 0 : CES-Atacama-LAT-small_patch-0-0 : Observing the scan-synchronous signal
Sigma is 76.438962 arcmin (0.022235 rad) 
-> fwhm is 180.000000 arcmin
Sigma is 0.000000 arcmin (0.000000 rad) 
-> fwhm is 0.000000 arcmin