HIDE: Hydrogen (HI) Data Emulator¶

HIDE is a package for simulating of a single dish radio telescope survey. As such, it takes healpix maps as inputs and processes them into TOD. The design is flexible and can be customized to different instruments and survey designs.

The HIDE package has been developed at ETH Zurich in the Software Lab of the Cosmology Research Group of the ETH Institute of Astronomy.
The development is coordinated on GitHub and contributions are welcome. The documentation of HIDE is available at readthedocs.org .
Contents:¶
Installation¶
The project is hosted on GitHub. Get a copy by running:
$ git clone https://github.com/cosmo-ethz/hide.git
Install the package like this:
$ cd hide
$ pip install -r requirements.txt
$ python setup.py install --user
Alternatively, if you want to develop new features:
$ cd hide
$ python setup.py develop --user
Usage¶
To use Hydrogen (HI) Data Emulator in a project execute the following on the command line:
$ hide --strategy-start=2016-03-21-00:00:00 --strategy-end=2016-03-21-23:59:00 --verbose=True hide.config.bleien7m
This will simulate one day of time-ordered-data from the Bleien 7m radio telescope.
To visualize 15 minutes of the generated data run this:
import matplotlib.pyplot as plt
import matplotlib
import h5py
with h5py.File("./2016/03/21/TEST_MP_PXX_20160321_000000.h5", "r") as fp:
tod = fp["P/Phase1"].value
time = fp["TIME"].value
plt.imshow(tod, aspect="auto",
extent=(time[0], time[-1],990, 1260),
cmap="gist_earth", norm=matplotlib.colors.LogNorm())
plt.colorbar()

hide Package¶
hide
Package¶
Subpackages¶
astro Package¶
gsm_point_src
Module¶
Created on Apr 22, 2016
In large parts a copy of astro_calibration_sources in seek by cchang.
Models taken from: Baars 1997, Hafez 2008, Benz 2009 All numbers divided by 2 to account for polarization.
Coordinates from wikipedia
author: seehars, jakeret
-
class
hide.astro.gsm_point_src.
AstroSource
(model, ra, dec)¶ Bases:
tuple
-
dec
¶ Alias for field number 2
-
model
¶ Alias for field number 0
-
ra
¶ Alias for field number 1
-
beam Package¶
beam
Package¶
-
hide.beam.
AxisSpec
¶ Beam definition
alias of
axis
airy
Module¶
Created on Apr 25, 2016
author: jakeret
gaussian
Module¶
Created on Dec 8, 2014
author: jakeret
gaussian_interp
Module¶
Created on Dec 8, 2014
author: jakeret
config Package¶
plugins Package¶
add_rfi
Module¶
Created on Dec 8, 2014
author: jakeret
add_rfi_phaseswitch
Module¶
Created on Feb 17, 2016
author: seehars
-
class
hide.plugins.add_rfi_phaseswitch.
Plugin
(ctx, **kwargs)[source]¶ Bases:
ivy.plugin.base_plugin.BasePlugin
Adds RFI to the time ordered data (phase switch).
-
hide.plugins.add_rfi_phaseswitch.
calcRFI
(background, amplitude, fraction, deltat, deltaf, exponent, enhance, nf, nt)[source]¶ Get time-frequency plane of RFI.
Parameters: - background – background level of data per channel
- amplitude – maximal amplitude of RFI per channel
- fraction – fraction of RFI dominated pixels per channel
- deltat – time scale of rfi decay (in units of pixels)
- deltaf – frequency scale of rfi decay (in units of pixels)
- exponent – exponent of rfi model (either 1 or 2)
- enhance – enhancement factor relative to fraction
- nf – number of frequency channels
- nt – number of time steps
Returns RFI: time-frequency plane of RFI
-
hide.plugins.add_rfi_phaseswitch.
getRFI
(background, amplitude, fraction, deltat, deltaf, exponent, enhance, frequencies, time, rfiday, damping)[source]¶ Get time-frequency plane of RFI.
Parameters: - background – background level of data per channel
- amplitude – maximal amplitude of RFI per channel
- fraction – fraction of RFI dominated pixels per channel
- deltat – time scale of rfi decay (in units of pixels)
- deltaf – frequency scale of rfi decay (in units of pixels)
- exponent – exponent of rfi model (either 1 or 2)
- enhance – enhancement factor relative to fraction
- frequencies – frequencies of tod in MHz
- time – time of day in hours of tod
- rfiday – tuple of start and end of RFI day
- damping – damping factor for RFI fraction during the RFI night
Returns RFI: time-frequency plane of RFI
-
hide.plugins.add_rfi_phaseswitch.
kernel
(deltaf, deltat, nf, nt, N, exponent)[source]¶ Convolution kernel for FFT convolution
Parameters: - deltaf – spread of RFI model in frequency
- deltat – spread of RFI model in time
- nf – number of frequencies
- nt – number of time steps
- N – size of kernel relative to deltaf, deltat
- exponent – exponent of RFI model (see logmodel)
Returns kernel: convolution kernel
-
hide.plugins.add_rfi_phaseswitch.
logmodel
(x, dx, exponent)[source]¶ - Model for the log of the RFI profile:
- -abs(x)/dx for exponent 1
- -(x/dx)^2 for exponent 2
Parameters: - x – grid on which to evaluate the profile
- dx – width of exponential
- exponent – exponent of (x/dx), either 1 or 2
Returns logmodel: log of RFI profile
background_noise
Module¶
Created on Dec 8, 2014
author: jakeret
coord_transform
Module¶
Created on Dec 8, 2014
author: jakeret
write_calibration
Module¶
Created on Mar 29, 2016
author: seehars
-
class
hide.plugins.write_calibration.
Plugin
(ctx, **kwargs)[source]¶ Bases:
ivy.plugin.base_plugin.BasePlugin
Writes the sources of the calibration days to disk
write_rfi
Module¶
Created on Feb 25, 2016
author: seehars
write_tod_phaseswitch
Module¶
Created on Dec 16, 2015
author: seehars
-
class
hide.plugins.write_tod_phaseswitch.
Plugin
(ctx, **kwargs)[source]¶ Bases:
ivy.plugin.base_plugin.BasePlugin
Writes the time ordered phase switch data to the file system
-
hide.plugins.write_tod_phaseswitch.
add_dataset
(grp, name, data)[source]¶ Adds a dataset to the group applying moderate compression
Parameters: - grp – The group
- name – Name of the dataset
- data – the actual data to be added
spectrometer Package¶
M9703A
Module¶
Created on Nov 9, 2015
author: jakeret
strategy Package¶
strategy
Package¶
drift_scan
Module¶
Created on Jan 15, 2015
author: jakeret
-
hide.strategy.drift_scan.
load_strategy
(ctx)[source]¶ Creates a scanning strategy that uses drift mode i.e. the telescope stares at the same position for 24 hours and then changes the altiude by a certain angle
Parameters: ctx – The ctx instance with the paramterization Returns strategy: A list of CoordSpec with the scanning strategy
scheduler
Module¶
Created on Mar 24, 2016
author: seehars
-
hide.strategy.scheduler.
load_strategy
(ctx)[source]¶ Creates a scanning strategy from a scheduler file.
Parameters: ctx – The ctx instance with the path to the scheduler file Returns strategy: A list of CoordSpec with the scanning strategy
-
hide.strategy.scheduler.
parse_schedule
(path, strategy_start)[source]¶ Parses a scheduler file :param path: the path to the scheduler file :param strategy_start: start date of the strategy
Returns schedule_entries: list of ScheduleEntry
-
hide.strategy.scheduler.
process_schedule
(schedule, step_size, strategy_start, strategy_end, obs)[source]¶ Processes a list of schedule entries :param schedule: the list of schedule entries :param step_size: the step size to use :param strategy_start: start date of the strategy :param strategy_end: end date of the strategy :param obs: telescope position
Returns strategy, calibration_days: a list of CoordSpec and a dict for the calibration days
scheduler_virtual
Module¶
Created on May 4, 2016
author: jakeret
utils Package¶
quaternion
Module¶
-
class
hide.utils.quaternion.
Rotator
(q)[source]¶ Bases:
object
Quaternion based rotator implementation for theta, phi
-
class
hide.utils.quaternion.
VecRotator
(q)[source]¶ Bases:
object
Quaternion based rotator implementation for theta, phi
-
hide.utils.quaternion.
mult
(p, q)[source]¶ Multiply arrays of quaternions, ndarray objects with 4 columns defined as x y z w see: http://en.wikipedia.org/wiki/Quaternions#Quaternions_and_the_geometry_of_R3
signal
Module¶
Created on Nov 10, 2015
author: jakeret
-
hide.utils.signal.
noisegen
(beta=0, N=8192)[source]¶ Noise will be generated that has spectral densities that vary as powers of inverse frequency, more precisely, the power spectra P(f) is proportional to 1 / fbeta for beta >= 0. When beta is 0 the noise is referred to white noise, when it is 2 it is referred to as Brownian noise, and when it is 1 it normally referred to simply as 1/f noise which occurs very often in processes found in nature.
The basic method involves creating frequency components which have a magnitude that is generated from a Gaussian white process and scaled by the appropriate power of f. The phase is uniformly distributed on 0, 2pi.
from http://paulbourke.net/fractals/noise/
Parameters: - beta –
- N – number of samples (can also be shape of array)
Returns out: the sampled noise
sphere
Module¶
Created on Dec 22, 2014
author: jakeret
-
class
hide.utils.sphere.
ArcKDTree
(theta, phi)[source]¶ Bases:
object
Wraps the scipy.spatial.cKDTree such that the tree can be used with spherical coords
-
query
(theta, phi, k=1, eps=0, p=2, distance_upper_bound=inf)[source]¶ Query the kd-tree for nearest neighbors using theta, phi :param theta: :param phi: :param k: :param eps: :param p: :param distance_upper_bound:
Returns d, i: The distances to the nearest neighbors, the locations of the neighbors in self.data.
-
-
hide.utils.sphere.
rotate_map
(Map, rotator, mask=None)[source]¶ Map is map in system A rotator is rotator from system B to A mask is a mask in system B returns new map in system B
-
hide.utils.sphere.
separation
(d1, a1, d2, a2)[source]¶ great circle distance http://en.wikipedia.org/wiki/Great-circle_distance#Computational_formulas
Parameters: - d1 – dec 1
- a1 – ra 1
- d2 – dec 2
:param a2:ra 2
Contributing¶
Contributions are welcome, and they are greatly appreciated! Every little bit helps, and credit will always be given.
You can contribute in many ways:
Types of Contributions¶
Report Bugs¶
If you are reporting a bug, please include:
- Your operating system name and version.
- Any details about your local setup that might be helpful in troubleshooting.
- Detailed steps to reproduce the bug.
Fix Bugs¶
Implement Features¶
Write Documentation¶
Hydrogen (HI) Data Emulator could always use more documentation, whether as part of the official Hydrogen (HI) Data Emulator docs, in docstrings, or even on the web in blog posts, articles, and such.
Submit Feedback¶
If you are proposing a feature:
- Explain in detail how it would work.
- Keep the scope as narrow as possible, to make it easier to implement.
- Remember that this is a volunteer-driven project, and that contributions are welcome :)
Pull Request Guidelines¶
Before you submit a pull request, check that it meets these guidelines:
- The pull request should include tests.
- If the pull request adds functionality, the docs should be updated. Put your new functionality into a function with a docstring, and add the feature to the list in README.rst.
- The pull request should work for Python 2.6, 2.7, and 3.3, and for PyPy. make sure that the tests pass for all supported Python versions.
Credits¶
Development Lead¶
- Joel Akeret <jakeret@phys.ethz.ch>
- Sebastian Seehars <seehars@phys.ethz.ch>
- Chihway Chang <chihway.chang@phys.ethz.ch>
Contributors¶
None yet. Why not be the first?
Feedback¶
If you have any suggestions or questions about Hydrogen (HI) Data Emulator feel free to email me at jakeret@phys.ethz.ch.
If you encounter any errors or problems with Hydrogen (HI) Data Emulator, please let me know!