API Reference

batch_tools

Tools for submitting to HPC systems: PBS or SLURM

class xfaster.batch_tools.JobArgumentParser(name=None, mem=None, time=None, workdir=None, outkey=None, **kwargs)[source]

Standardized way to add job submission arguments to a script.

Keyword arguments are used to fix values for parameters not needed by a script. The corresponding command line arguments will not be added. See the .opt_list attribute for a complete list.

Parameters:
  • name (string) – This will be the name of the jobs.

  • mem (float) – Memory per job in GB. Will scale up for grouped jobs.

  • time (float) – Time per job in hours.

  • workdir (string) – Fixed location to use as working directory (ie place for job scripts and logs).

  • outkey (string) – Alternative to workdir. The key/name of an argument added to normal argparse that indicates output path. A “logs” subfolder of that path will be used as workdir.

Example

>>> jp = sa.batch.JobArgumentParser("some_serial_job", mem=4, time=1.5,
            omp_threads=1, workdir="/path/to/logs")
>>> AP = argparse.ArgumentParser(description="Do some job")
>>> AP.add_argument(...)  # other non-job arguments
>>> jp.add_arguments(AP)
>>> args = AP.parse_args()
>>> jp.set_job_opts(args)
>>> jobs = [...]  # list of commands to run in jobs
>>> jp.submit(jobs)
add_arguments(parser=None, add_group=True)[source]

Add job submission arguments to an argparse.ArgumentParser.

Parameters:
  • parser (argparse.ArgumentParser) – The parser to which to add arguments. If None, a new parser will be made and returned.

  • add_group (bool or string) – Whether to add job submit options in an argument group. If a string, use as description for the group.

Returns:

parser – The updated argument parser object.

Return type:

argparse.ArgumentParser

pop_job_opts(args_dict, pop_submit=True)[source]

Pop all of the job-related options from a dictionary of arguments

Parameters:
  • args_dict (dict) – The dictionary to pop from

  • pop_submit (bool) – Whether to also pop an argument named “submit”

set_job_opts(args, load_defaults=True, **kwargs)[source]

Set job submission options based on parsed arguments.

Keyword arguments will be passed to update. Can be used to override particular job submission options. Any argument to the batch_sub or batch_group functions can be overridden in this way.

Parameters:
  • args (argparse.Namespace or dict) – The parsed command line arguments (from argparse.ArgumentParser.parse_args()).

  • load_defaults (bool) – Whether to automatically load the default value for options

submit(jobs, **kwargs)[source]

Submit jobs based on the parsed arguments. Must be called after set_job_opts.

Keyword arguments will be passed to update. Can be used to override particular job submission options. Any argument to the batch_sub or batch_group functions can be overridden in this way.

Parameters:

jobs (string or list of strings) – The job(s) to submit.

Return type:

List of job IDs for submitted jobs.

update(**kwargs)[source]

Update particular job submission options.

Keyword arguments can be used to override particular job submission options. Any argument to the batch_sub or batch_group functions can be overridden in this way.

xfaster.batch_tools.batch_group(cmds, group_by=1, serial=False, *args, **kwargs)[source]

Create and submit SLURM job scripts for a group of similar commands. The commands can be grouped together into larger single jobs that run them in parallel on multiple processors on a node.

Keyword arguments are passed on to the batch_sub function. These will be applied to EACH job. For example, using nodes="1:ppn=8" with group_by=8 and 16 elements in cmds will result in 2 jobs, each using 8 processors on 1 node.

Parameters:
  • cmds (list) – The commands to run. The commands themselves should be a string, or a list of tokens, as per the batch_sub function.

  • group_by (int, optional) – The number of commands to group together into a single job. Does not balance well when len(cmds)%group_by != 0 Eg. on scinet use group_by=8 to efficiently use whole nodes.

  • serial (bool, optional) – Set to True to run cmds sequentially, rather than starting them all in parallel. This will also work with OpenMP parallel jobs.

  • args (arb) – Additional arguments passed to batch_sub

  • kwargs (arb) – Additional arguments passed to batch_sub

Returns:

jobids – List of job IDs used by the scheduler

Return type:

list of strings

xfaster.batch_tools.batch_sub(cmd, name=None, mem=None, nodes=None, node_list=None, ppn=None, cput=None, wallt=None, output=None, error=None, queue=None, dep_afterok=None, workdir=None, batch_args=[], omp_threads=None, env_script=None, env=None, nice=None, echo=True, delete=True, submit=True, scheduler='slurm', debug=False, exclude=None, verbose=False, srun=None, srun_args='')[source]

Create and submit a SLURM job.

Parameters:
  • cmd (string or list of strings) – A command sequence to run via SLURM. The command will be inserted into a qsub submission script with all of the options specified in the remaining arguments.

  • name (string, optional) – Name of the job.

  • mem (float or string, optional) – Amount of memory to request for the job. float values in GB. Or pass a string (eg ‘4gb’) to use directly.

  • nodes (int or string, optional) – Number of nodes to use in job If using SLURM and a string, will overwrite node_list if None

  • node_list (string or list of strings) – List of nodes that can be used for job. SLURM-only.

  • ppn (int, optional) – Number of processes per node

  • cput (string or float or datetime.timedelta, optional) – Amount of CPU time requested. String values should be in the format HH:MM:SS, e.g. ‘10:00:00’. Numerical values are interpreted as a number of hours.

  • wallt (string or float or datetime.timedelta, optional) – Amount of wall clock time requested. String values should be in the format HH:MM:SS, e.g. ‘10:00:00’. Numerical values are interpreted as a number of hours.

  • output (string, optional) – Scheduler standard output filename.

  • error (string, optional) – Scheduler error output filename.

  • queue (string, optional) – The name of the queue to which to submit jobs

  • dep_afterok (string or list of strings) – Dependency. Job ID (or IDs) on which to wait for successful completion, before starting this job

  • workdir (string, optional) – Directory from where the script will be submitted. This is where the output and error files will be created by default. Default: current directory.

  • batch_args (string or list of strings, optional) – Any additional arguments to pass to slurm.

  • omp_threads (int, optional) – Number of OpenMP threads to use per process

  • env_script (string, optional) – Path to script to source during job script preamble For loading modules, setting environment variables, etc

  • env (dict, optional) – Dictionary of environment variables to set in job script

  • nice (int, optional) – Adjust scheduling priority (SLURM only). Unprivileged users can only increase niceness (lower job priority).

  • echo (bool, optional) – Whether to use bash “set -x” in job script to echo commands to stdout.

  • delete (bool, optional) – If True, delete the submit script upon job submission.

  • submit (bool, optional) – If True (default) submit the job script once create. Will override the default option when False, to keep the script

  • scheduler (string, optional) – Which scheduler system to write a script for. Only “slurm” is currently supported.

  • debug (bool, optional) – If True, print the contents of the job script to stdout for debugging.

  • exclude (string or list of strings) – List of nodes that will be excluded for job. SLURM-only.

  • verbose (bool, optional) – Print the working directory, and the job ID if submitted successfully.

  • srun (bool, optional) – Invoke jobs with srun. Defaults to True when using SLURM, set to False to disable. Use srun_args to pass extra arguments to this command.

  • srun_args (string, optional) – If srun=True, additional arguments to pass to the srun command.

Returns:

jobid – The ID of the submitted job.

Return type:

string

Example

>>> jobid = batch_sub("echo Hello", name="testing", nodes="1:ppn=1",
... cput='1:00:00', mem='1gb')
>>> print(jobid)
221114.feynman.princeton.edu
>>> print(open('testing.o221114','r').read())
Hello
xfaster.batch_tools.get_job_id()[source]

Get the current job ID, based on job environment.

Returns:

jobid – Job ID.

Return type:

str

xfaster.batch_tools.get_job_logfile()[source]

Generate a path to use for the output log, based on job environment

Returns:

logfile – Path to log file.

Return type:

str

parse_tools

Module for parsing intermediate and output XFaster files.

Tools for parsing XFasters outputs and intermediate data structures

xfaster.parse_tools.arr_to_dict(arr, ref_dict)[source]

Transform an array of data into a dictionary keyed by the same keys in ref_dict, with data divided into chunks of the same length as in ref_dict. Requires that the length of the array is the sum of the lengths of the arrays in each entry of ref_dict. The other dimensions of the input array and reference dict can differ.

Parameters:
  • arr (array) – Input array to be transformed into dictionary.

  • ref_dict (dict) – Reference dictionary containing the keys used to construct the output dictionary.

Returns:

out – Dictionary of values from arr keyed with keys from ref_dict.

Return type:

dict

xfaster.parse_tools.bin_spec_simple(qb, cls_shape, bin_def, inv_fish=None, lfac=True)[source]

Compute binned output spectra and covariances by averaging the shape spectrum over each bin, and applying the appropriate qb bandpower amplitude. NB: this does not use the true window functions to compute bandpowers, and the results should be treated as an approximation.

Parameters:
  • qb (dict) – Bandpower amplitudes for each spectrum bin.

  • cls_shape (dict) – Shape spectrum

  • bin_def (dict) – Bin definition dictionary

  • inv_fish (array_like, (nbins, nbins)) – Inverse fisher matrix for computing the bin errors and covariance. If not supplied, these are not computed.

  • lfac (bool) – If False, return binned C_l spectrum rather than the default D_l

Returns:

  • cb (dict of arrays) – Binned spectrum

  • dcb (dict of arrays) – Binned spectrum error, if inv_fish is not None

  • ellb (dict of arrays) – Average bin center

  • cov (array_like, (nbins, nbins)) – Binned spectrum covariance, if inv_fish is not None

  • qb2cb (dict) – The conversion factor from qb to cb, computed by averaging over the input shape spectrum.

xfaster.parse_tools.corr_index(idx, n)[source]

Gets the index of the auto spectrum when getting all pairwise combinations of n maps.

Parameters:
  • idx (int) – The index of the map in the list of maps being looped through.

  • n (int) – The number of maps being looped through.

Returns:

index – Index of auto spectrum

Return type:

int

xfaster.parse_tools.dict_to_arr(d, out=None, flatten=False)[source]

Transform ordered dict into an array, if all items are same shape

If not all items are the same shape, eg, for qb, or if flatten=True, flatten everything into a vector

Parameters:
  • d (dict) – Dictionary to transform into an array.

  • out (array) – If not None, the starting array onto which to stack the arrays contained in dictionary d.

  • flatten (bool) – If True, return flattened vector instead of array.

Returns:

out – The array containing stacked elements of the arrays contained in dictionary d.

Return type:

array

xfaster.parse_tools.dict_to_dmat(dmat_dict, pol=None)[source]

Take a dmat dictionary and return the right shaped Dmat matrix: (Nmaps * 3, Nmaps * 3, lmax + 1) if pol else (Nmaps, Nmaps, lmax + 1)

Parameters:

dmat_dict (dict) – Dictionary containing the model covariance terms.

Returns:

Dmat – Dmat total model covariance matrix.

Return type:

arr

xfaster.parse_tools.dict_to_dsdqb_mat(dsdqb_dict, bin_def)[source]

Take a dSdqb dictionary and return the right shaped dSdqb matrix: (Nmaps * 3, Nmaps * 3, nbins_cmb+nbins_fg+nbins_res, lmax + 1) if pol else first two dimensions are Nmaps.

Parameters:
  • dsdqb_dict (dict) – Dictionary containing the terms for the derivative of the signal model, S, w.r.t. the qb parameters.

  • bin_def (dict) – Dictionary containing the bin edges for each qb value fit.

Returns:

dsdqb_mat – Signal derivative matrix in the expected shape for matrix multiplication in the Fisher iteration.

Return type:

arr

xfaster.parse_tools.dict_to_index(d)[source]

Construct a dictionary of (start, stop) indices that correspond to the location of each sub-array when the dict is converted to a single array using dict_to_arr.

Parameters:

d (dict) – Input dictionary.

Returns:

index – Dictionary containing location of sub-arrays corresponding to keys.

Return type:

dict

Examples

To use this function to index into a (nbins, nbins) array, create the index dictionary:

>>> bin_def = OrderedDict((k, np.array([[2, 27], [27, 52]]))
...                       for k in ['cmb_tt', 'cmb_ee', 'cmb_bb'])
>>> bin_index = dict_to_index(bin_def)
>>> bin_index
OrderedDict([('cmb_tt', (0, 2)),
             ('cmb_ee', (2, 4)),
             ('cmb_bb', (4, 6))])

To extract the TT bins from the fisher matrix:

>>> fisher = np.random.randn(12, 12)
>>> sl_tt = slice(*bin_index['cmb_tt'])
>>> fisher_tt = fisher[sl_tt, sl_tt]

To extract all the CMB bins from the fisher matrix:

>>> sl_cmb = slice(bin_index['cmb_tt'][0], bin_index['cmb_bb'][1])
>>> fisher_cmb = fisher[sl_cmb, sl_cmb]
xfaster.parse_tools.expand_qb(qb, bin_def, lmax=None)[source]

Expand a qb-type array to an ell-by-ell spectrum using bin_def.

Parameters:
  • qb (array_like, (nbins,)) – Array of bandpower deviations

  • bin_def (array_like, (nbins, 2)) – Array of bin edges for each bin

  • lmax (int, optional) – If supplied, limit the output spectrum to this value. Otherwise the output spectrum extends to include the last bin.

Returns:

cl – Array of expanded bandpowers

Return type:

array_like, (lmax + 1,)

xfaster.parse_tools.fix_data_roots(data, mode='save', root=None, root2=None)[source]

Remove or apply the data root to a set of file paths in an output checkpoint file. Operations are performed in-place on the input dictionary.

Parameters:
  • data (dict) – Dictionary of file data to update, including lists of map files and file roots, as returned by XFaster.get_files and related functions.

  • mode (str) – If 'save', remove the appropriate data root from any file paths, so that only relative paths are stored to disk. If 'load', add the appropriate data root to any file paths, so that all paths point to existing files within the given data roots.

  • root (str) – The data root (and second data root, for null tests) to be removed or applied to the file path items in data.

  • root2 (str) – The data root (and second data root, for null tests) to be removed or applied to the file path items in data.

Returns:

data – Updated data dictionary.

Return type:

dict

xfaster.parse_tools.load_and_parse(filename, check_version=True)[source]

Load a .npz data file from disk and parse all the fields it contains. Includes handling of backward compatibility to older file versions on disk.

Returns a dictionary of parsed fields.

Parameters:
  • filename (str) – Path to numpy data file on disk.

  • check_version (bool) – If True, check the data file version and apply any necessary updates to the latest version.

Returns:

data – Data dictionary loaded from disk.

Return type:

dict

xfaster.parse_tools.load_compat(*args, **kwargs)[source]

Load and decode a numpy archive file from disk.

Backward compatible with python2 data files.

Parameters:
  • args (key/value pairs) – Passed to np.load.

  • kwargs (key/value pairs) – Passed to np.load.

Returns:

out – Dictionary of info from numpy archive file

Return type:

dict

xfaster.parse_tools.load_pickle_compat(filename)[source]

Load a pickle file from the given filename. Ensure that the file is open in mode ‘rb’ (required for python3), and that the encoding is set to ‘latin1’ in python3.

Parameters:

filename (str) – Path to pickled output file to read.

Returns:

f – Unpickled file.

Return type:

dict

xfaster.parse_tools.num_corr(n)[source]

Returns how many cross spectra there are if there are n total maps.

Parameters:

n (int) – Number of maps.

Returns:

nxspec – Number of cross spectra.

Return type:

int

xfaster.parse_tools.num_maps(n)[source]

Returns how many maps there are if there are n total cross spectra.

Parameters:

n (int) – Number of cross spectra.

Returns:

nmaps – Number of maps.

Return type:

int

xfaster.parse_tools.save(output_file, **data)[source]

Save a numpy archive file to disk.

Parameters:
  • filename (str) – Path to output npz file.

  • data (dict) – Dictionary of data to store.

xfaster.parse_tools.spec_index(spec=None)[source]

Return the matrix indices of the given spectrum within a 3x3 matrix. If spec is None, return a dictionary of such indices keyed by spectrum.

Parameters:

spec (str) – Which spectrum to return index for. If None, return dict of all.

Returns:

index – Dictionary of indices if spec in None, or list of indices if spec is provided.

Return type:

dict or list

xfaster.parse_tools.spec_mask(spec=None, nmaps=1)[source]

Return a mask for extracting spectrum terms from a matrix of shape (3 * nmaps, 3 * nmaps). If spec is None, returns a dictionary of masks keyed by spectrum.

Parameters:
  • spec (str) – Which spectrum to return mask for. If None, return dict of all masks.

  • nmaps (int) – Number of maps used for the cross-spectrum analysis.

Returns:

spec_mask – Dictionary of masks if spec in None, or (3 * nmaps, 3 * nmaps) array that is 1 in elements corresponding to spec if spec is provided.

Return type:

dict or arr

xfaster.parse_tools.tag_pairs(tags, index=False)[source]

Return an OrderedDict whose keys are pairs of tags in the format “tag1:tag2” and whose values are a tuple of the two tags used to construct each key, or a tuple of the indices of the two tags in the original tag list, if index is True. If index is a list, then it should be a list the same length as tags, and the tuple is populated by indexing into index using the two indices of the tags in the original tag list.

Parameters:
  • tags (list of strings) – Map tags from which to construct cross-spectrum keys like “tag1:tag2”.

  • index (bool) – If True, make values in dictionary the indices of the map tags, rather than the tags themselves.

Returns:

pairs – Dictionary whose keys are pairs of tags in the format “tag1:tag2” and whose values are a tuple of the two tags used to construct the key, or their indices, if index=True.

Return type:

OrderedDict

Example

>>> tags = ['a', 'b']
>>> tag_pairs(tags)
OrderedDict([('a:a', ('a', 'a')), ('a:b', ('a', 'b')), ('b:b', ('b', 'b'))])
>>> tag_pairs(tags, index=True)
OrderedDict([('a:a', (0, 0)), ('a:b', (0, 1)), ('b:b', (1, 1))])
>>> tag_pairs(tags, index=['c', 'd'])
OrderedDict([('a:a', ('c', 'c')), ('a:b', ('c', 'd')), ('b:b', ('d', 'd'))])
xfaster.parse_tools.unique_tags(tags)[source]

If map tags are repeated (eg, two 150 maps in different chunk subdirectories), return a list modifying them with an index

Parameters:

tags (list of strings) – List of original map tags.

Returns:

new_tags – List of map tags where repeated tags are modified to be <tag>_<index>, with unique indices.

Return type:

list of strings

xfaster_class

The main class where the algorithmic steps are defined

class xfaster.xfaster_class.XFaster(config, output_root=None, output_tag=None, verbose='notice', debug=False, checkpoint=None, alm_pixel_weights=False, alm_iter=None, add_log=False)[source]

Initialize an XFaster instance for computing binned power spectra using a set of data maps along with signal and noise simulations.

Parameters:
  • config (string) – Configuration file. If path doesn’t exist, assumed to be in xfaster/config/<config>

  • output_root (string) – Path to data output directory

  • output_tag (string) – Tag to use for output data. Results are typically stored in the form <output_root>/<output_tag>/<name>_<output_tag>.npz

  • verbose (string) – Verbosity level to use for log messages. Can be one of [“critical”, “error”, “warning”, “notice”, “info”, “debug”, “all”].

  • debug (bool) – Store extra data in output files for debugging.

  • checkpoint (string) – If output data from this step forward exist on disk, they are are re-computed rather than loading from file. Options are [‘files’, ‘masks’, ‘kernels’, ‘sims_transfer’, ‘shape_transfer’, ‘transfer’, ‘sims’, ‘beams’, ‘data’, ‘sim_data’, ‘template_noise’, ‘shape’, ‘bandpowers’, ‘beam_errors’, ‘likelihood’].

  • alm_pixel_weights (bool) – If True, set the use_pixel_weights option to True when computing map Alms using healpy.map2alm. If False, sets the use_weights option to True instead. Note: pixel weights ensure accurate Alm computation, but can only be used for analyses where lmax < 1.5 * nside.

  • alm_iter (int) – If given, set the iter option to the given value when computing map Alms using healpy.map2alm. Using more iterations improves the accuracy of the output Alms. If not set, uses the default number of iterations (3) as defined in healpy.

  • add_log (bool) – If True, write log output to a file instead of to STDOUT. The log will be in <output_root>/run_<output_tag>.log.

alm2cl(m1, m2=None, lmin=2, lmax=None, symmetric=True)[source]

Wrapper for healpy.alm2cl.

Parameters:
  • m1 (array_like) – Masked alms for map1

  • m2 (array_like) – Masked alms for map2

  • lmin (int, default: 2) – The minimum ell bin to include in the output Cls. All ell bins below this are nulled out.

  • lmax (int) – The maximum ell bin to compute. If None, this is set to the lmax value with which the class was initialized.

  • symmetric (bool, default: True) – If True, the average cross spectrum of (m1-x-m2 + m2-x-m1) / 2. is computed.

Returns:

cls – Cross-spectrum of m1-x-m2.

Return type:

array_like

apply_mask(m, mask)[source]

Apply the input mask to the data map, in place.

If the map is polarized, the appropriate mask is applied to the polarization data, depending on whether the mask is also polarized.

Parameters:
  • m (array_like) – Input map (T/Q/U if polarized, T-only if not) This array is modified in place.

  • mask (array_like) – Mask to apply (T/P if polarized, T-only if not)

Returns:

m – Input map multiplied by mask.

Return type:

array_like

bin_cl_template(cls_shape=None, map_tag=None, transfer=False, beam_error=False, use_precalc=True)[source]

Compute the Cbl matrix from the input shape spectrum.

This method requires beam windows, kernels and transfer functions (if transfer is False) to have been precomputed.

Parameters:
  • cls_shape (array_like) – The shape spectrum to use. This can be computed using get_signal_shape or otherwise.

  • map_tag (str) – If supplied, the Cbl is computed only for the given map tag (or cross if map_tag is map_tag1:map_tag2). Otherwise, it is computed for all maps and crosses.

  • transfer (bool) – If True, this assumes a unity transfer function for all bins, and the output Cbl is used to compute the transfer functions that are then loaded when this method is called with transfer = False.

  • beam_error (bool) – If True, use beam error envelope instead of beam to get cbls that are 1 sigma beam error envelope offset of signal terms.

  • use_precalc (bool) – If True, load pre-calculated terms stored from a previous iteration, and store for a future iteration. Otherwise, all calculations are repeated.

Returns:

cbl – The Cbl matrix, indexed by component and spectrum, then by map cross, e.g. cbl['cmb_tt']['map1:map2']. E/B mixing terms are stored in elements cbl['cmb_ee_mix'] and cbl['cmb_bb_mix'], and unmixed terms are stored in elements cbl['cmb_ee'] and cbl['cmb_bb'].

Return type:

dict of arrays (num_bins, 2, lmax + 1)

checkpoint_tree = {'bandpowers': ['likelihood'], 'beam_errors': ['likelihood'], 'beams': ['transfer'], 'data': ['bandpowers'], 'files': ['masks'], 'kernels': ['transfer'], 'masks': ['kernels', 'sims_transfer', 'sims', 'data', 'sim_data', 'template_noise'], 'shape': ['bandpowers'], 'shape_transfer': ['transfer'], 'sim_data': ['bandpowers'], 'sims': ['bandpowers'], 'sims_transfer': ['transfer'], 'template_noise': ['bandpowers'], 'transfer': ['bandpowers']}
checkpoints = ['files', 'masks', 'kernels', 'sims_transfer', 'shape_transfer', 'transfer', 'sims', 'beams', 'data', 'sim_data', 'template_noise', 'shape', 'bandpowers', 'beam_errors', 'likelihood']
clear_precalc()[source]

Clear variables pre-computed with fisher_precalc.

data_version = 4
do_qb2cb(qb, inv_fish, wbl)[source]

Compute binned output spectra and covariances by averaging the shape spectrum over each bin, and applying the appropriate qb bandpower amplitude.

This method is used internally by fisher_calc, and requires bin definitions to have been pre-loaded using get_bin_def or get_transfer.

Parameters:
  • qb (dict) – Bandpower amplitudes for each spectrum bin.

  • inv_fish (array_like, (nbins, nbins)) – Inverse fisher matrix for computing the bin errors and covariance.

  • wbl (dict) – Window functions for each qb

Returns:

  • cb (dict of arrays) – Binned spectrum

  • dcb (dict of arrays) – Binned spectrum error

  • ellb (dict of arrays) – Average bin center

  • cov (array_like, (nbins, nbins)) – Binned spectrum covariance

  • qb2cb (dict) – The conversion matrix from qb to cb for each spectrum component, computed from the qb window functions

  • wbl_cb (dict) – Window functions for each cb

fisher_calc(qb, cbl, cls_input, cls_noise=None, cls_debias=None, cls_model=None, cond_noise=None, cond_criteria=None, likelihood=False, like_lmin=2, like_lmax=None, delta_beta_prior=None, null_first_cmb=False, use_precalc=True, windows=False, inv_fish=None)[source]

Re-compute the Fisher matrix and qb amplitudes based on input data spectra. This method is called iteratively by fisher_iterate, and requires bin definitions precomputed by get_bin_def or get_transfer.

Parameters:
  • qb (OrderedDict) – Bandpower amplitudes, typically computed in a previous call to this method.

  • cbl (OrderedDict) – Cbl matrix computed by bin_cl_template for a given shape spectrum.

  • cls_input (OrderedDict) – Input spectra. If computing a transfer function, this is the average cls_signal. If computing a null test, this is cls_data_null, and otherwise it is cls_data, for a single map or several input maps.

  • cls_noise (OrderedDict) – If supplied, the noise spectrum is applied to the model spectrum.

  • cls_debias (OrderedDict) – If supplied, the debias spectrum is subtracted from the input.

  • cls_model (OrderedDict) – Unbinned model spectrum computed from cbl

  • cond_noise (float) – The level of regularizing noise to add to EE and BB diagonals.

  • cond_criteria (float) – The maximum condition number allowed for Dmat1 to be acceptable for taking its inverse.

  • likelihood (bool) – If True, return the likelihood for the given input bandpowers, shapes and data spectra. Otherwise, computes output bandpowers and the fisher covariance for a NR iteration.

  • like_lmin (int) – The minimum ell value to be included in the likelihood calculation

  • like_lmax (int) – The maximum ell value to be included in the likelihood calculation

  • delta_beta_prior (float) – The width of the prior on the additive change from beta_ref. If you don’t want the code to fit for a spectral index different from beta_ref, set this to be None.

  • null_first_cmb (bool) – Keep first CMB bandpowers fixed to input shape (qb=1).

  • use_precalc (bool) – If True, load pre-calculated terms stored from a previous iteration, and store for a future iteration. Otherwise, all calculations are repeated.

  • windows (bool) – If True, return W_bl window functions for each CMB qb.

  • inv_fish (array_like) – Inverse Fisher matrix. If provided, don’t need to recompute. Useful if just getting bandpower window functions.

Returns:

  • qb (OrderedDict) – New bandpower amplitudes

  • inv_fish (array_like) – Inverse Fisher correlation matrix over all bins

  • – or –

  • likelihood (scalar) – Likelihood of the given input parameters.

  • – or –

  • windows (OrderedDict) – qb window functions

fisher_iterate(cbl, map_tag=None, iter_max=200, converge_criteria=0.005, qb_start=None, transfer=False, save_iters=False, null_first_cmb=False, delta_beta_prior=None, cond_noise=None, cond_criteria=None, qb_only=False, like_profiles=False, like_profile_sigma=3.0, like_profile_points=100, file_tag=None)[source]

Iterate over the Fisher calculation to compute bandpower estimates assuming an input shape spectrum.

Parameters:
  • cbl (OrderedDict) – Cbl matrix computed from an input shape spectrum

  • map_tag (str) – If not None, then iteration is performed over the spectra corresponding to the given map, rather over all possible combinations of map-map cross-spectra. In this case, the first dimension of the input cbl must be of size 1 (this is done automatically by calling bin_cl_template(..., map_tag=<map_tag>).

  • iter_max (int) – Maximum number of iterations to perform. if this limit is reached, a warning is issued.

  • converge_criteria (float) – Maximum fractional change in qb that indicates convergence and stops iteration.

  • qb_start (OrderedDict) – Initial guess at qb bandpower amplitudes. If None, unity is assumed for all bins.

  • transfer (bool) – If True, the input Cls passed to fisher_calc are the average of the signal simulations, and noise cls are ignored. If False, the input Cls are either cls_data_null (for null tests) or cls_data. See get_masked_data for how these are computed. The input noise is similarly either cls_noise_null or cls_noise.

  • save_iters (bool) – If True, the output data from each Fisher iteration are stored in an individual npz file.

  • null_first_cmb (bool) – Keep first CMB bandpowers fixed to input shape (qb=1).

  • delta_beta_prior (float) – The width of the prior on the additive change from beta_ref. If you don’t want the code to fit for a spectral index different from beta_ref, set this to be a very small value (O(1e-10)).

  • cond_noise (float) – The level of regularizing noise to add to EE and BB diagonals.

  • cond_criteria (float) – The maximum condition number allowed for Dmat1 to be acceptable for taking its inverse.

  • qb_only (bool) – If True, compute just maximum likelihood qb’s, and skip calculation of window functions and Cb bandpowers.

  • like_profiles (bool) – If True, compute profile likelihoods for each qb, leaving all others fixed at their maximum likelihood values. Profiles are computed over a range +/–sigma as estimated from the diagonals of the inverse Fisher matrix.

  • like_profile_sigma (float) – Range in units of 1sigma over which to compute profile likelihoods

  • like_profile_points (int) – Number of points to sample along the likelihood profile

  • file_tag (string) – If supplied, appended to the bandpowers filename.

Returns:

data – The results of the Fisher iteration process. This dictionary contains the fields:

qb : converged bandpower amplitudes
cb : output binned spectrum
dcb : binned spectrum errors
ellb : bin centers
cov : spectrum covariance
inv_fish : inverse fisher matrix
fqb : fractional change in qb for the last iteration
qb2cb : conversion factor
cbl : Cbl matrix
cls_model : unbinned model spectrum computed from Cbl
bin_def : bin definition array
cls_obs : observed input spectra
cls_noise : noise spectra
cls_shape : shape spectrum
iters : number of iterations completed

If transfer is False, this dictionary also contains:

qb_transfer : transfer function amplitudes
transfer : ell-by-ell transfer function
nbins_res : number of residual bins

Return type:

dict

Notes

This method stores outputs to files with name ‘transfer’ for transfer function runs (if transfer = True), otherwise with name ‘bandpowers’. Outputs from each individual iteration, containing only the quantities that change with each step, are stored in separate files with the same name (and different index).

fisher_precalc(cbl, cls_input, cls_debias=None, likelihood=False, windows=False)[source]

Pre-compute the D_ell and signal derivative matrices necessary for fisher_calc from the input data spectra. This method requires bin definitions precomputed by get_bin_def or get_transfer.

Parameters:
  • cbl (OrderedDict) – Cbl dict computed by bin_cl_template for a given shape spectrum.

  • cls_input (OrderedDict) – Input spectra. If computing a transfer function, this is the average cls_signal. If computing a null test, this is cls_data_null, and otherwise it is cls_data, for a single map or several input maps.

  • cls_debias (OrderedDict) – If supplied, the debias spectrum is subtracted from the input.

  • likelihood (bool) – If True, compute just Dmat_obs_b. Otherwise, Dmat_obs and dSdqb_mat1 are also computed.

  • windows (bool) – If True, compute dSdqb and Mll for constructing window functions.

Returns:

  • Dmat_obs (OrderedDict) – De-biased D_ell matrix from cls_input

  • Dmat_obs_b (OrderedDict) – Biased D_ell matrix from cls_input (for likelihood)

  • dSdqb_mat1 (OrderedDict) – Signal derivative matrix from Cbl

  • Mmat (OrderedDict) – Mode mixing matrix (Kll’ * Fl * Bl^2) for constructing window functions.

  • .. note:: the output arrays are also stored as attributes of the

  • parent object to avoid repeating the computation in fisher_calc

force_rerun_children(checkpoint)[source]

Trigger rerunning steps that depend on this checkpoint.

get_bandpowers(map_tag=None, converge_criteria=0.005, iter_max=200, return_qb=False, save_iters=False, delta_beta_prior=None, cond_noise=None, cond_criteria=None, null_first_cmb=False, return_cls=False, qb_only=False, like_profiles=False, like_profile_sigma=3.0, like_profile_points=100, file_tag=None)[source]

Compute the maximum likelihood bandpowers of the data, assuming a given input spectrum shape. Requires the transfer function to have been computed and loaded using get_transfer.

Parameters:
  • map_tag (string) – If not None, then iteration is performed over the spectra corresponding to the given map, rather over all possible combinations of map-map cross-spectra. In this case, the first dimension of the input cbl must be of size 1 (this is done automatically by calling bin_cl_template(..., map_tag=<map_tag>).

  • converge_criteria (float) – Maximum fractional change in qb that indicates convergence and stops iteration.

  • iter_max (int) – Maximum number of iterations to perform. if this limit is reached, a warning is issued.

  • return_qb (bool) – If True, only the maximum likelihood qb values are returned. Otherwise, the complete output dictionary is returned.

  • save_iters (bool) – If True, the output data from each Fisher iteration are stored in an individual npz file.

  • delta_beta_prior (float) – The width of the prior on the additive change from beta_ref. If you don’t want the code to fit for a spectral index different from beta_ref, set this to be None.

  • cond_noise (float) – The level of regularizing noise to add to EE and BB diagonals.

  • cond_criteria (float) – The maximum condition number allowed for Dmat1 to be acceptable for taking its inverse.

  • null_first_cmb (bool) – Keep first CMB bandpowers fixed to input shape (qb=1).

  • return_cls (bool) – If True, return C_ls rather than D_ls

  • qb_only (bool) – If True, compute just maximum likelihood qb’s, and skip calculation of window functions and Cb bandpowers.

  • like_profiles (bool) – If True, compute profile likelihoods for each qb, leaving all others fixed at their maximum likelihood values. Profiles are computed over a range +/–sigma as estimated from the diagonals of the inverse Fisher matrix.

  • like_profile_sigma (float) – Range in units of 1sigma over which to compute profile likelihoods

  • like_profile_points (int) – Number of points to sample along the likelihood profile

  • file_tag (string) – If supplied, appended to the bandpowers filename.

Returns:

  • data (dict) – Dictionary of maximum likelihood quantities, as output by fisher_iterate.

  • – or –

  • qb, inv_fish (array_like) – Maximum likelihood bandpower amplitudes and fisher covariance.

Notes

This method is called at the ‘bandpowers’ checkpoint, and loads or saves a data dictionary named ‘bandpowers’ with the quantities returned by fisher_iterate.

get_beam_errors()[source]

Get error envelope to multiply beam by (so, to get beam + 2 sigma error, do beam * (1 + 2 * beam_error))

get_beams(pixwin=True)[source]

Return beam window functions for all input map tags.

Parameters:

pixwin (bool) – If True, the pixel window function for the map nside is applied to the Gaussian beams.

Returns:

beam_windows – A dictionary of beam window function arrays (3 x lmax+1 if pol, 1 x lmax+1 if not) for each map tag

Return type:

dict

Notes

This method is called at the ‘beams’ checkpoint and loads or saves a dictionary containing just the beam_windows key to disk.

get_bin_def(lmin=2, lmax=500, pol=True, pol_mask=True, tbeb=False, bin_width=25, weighted_bins=False, residual_fit=False, res_specs=None, bin_width_res=25, foreground_fit=False, beta_fit=False, bin_width_fg=25, lmin_fg=None, lmax_fg=None)[source]

Construct the bin definition array that defines the bins for each output spectrum.

Parameters:
  • lmin (int) – Minimum ell for binned spectra

  • lmax (int) – The maximum multipole for which spectra are computed

  • pol (bool) – If True, polarized spectra are computed from the input maps

  • pol_mask (bool) – If True, two independent masks are applied to every map: one for T maps and one for Q/U maps.

  • tbeb (bool) – If True, EB and TB bins are constructed so that these spectra are computed by the estimator. Otherwise, these spectra are fixed at zero.

  • bin_width (int or list of ints) – Width of each spectrum bin. If a scalar, the same width is applied to all cross spectra. Otherwise, must be a list of up to six elements, listing bin widths for the spectra in the order (TT, EE, BB, TE, EB, TB).

  • weighted_bins (bool) – If True, use an lfac-weighted binning operator to construct Cbls. By default, a flat binning operator is used.

  • residual_fit (bool) – If True, fit for (compute bandpower amplitudes for) several wide bins of excess noise.

  • res_specs (list of strings) – List of spectra which are to be included in the residual fit. Can be individual spectra (‘tt’, ‘ee’, ‘bb’), or ‘eebb’ to fit for EE and BB residuals simultaneously. If not supplied, this defaults to [‘eebb’] for polarized maps, and [‘tt’] for unpolarized maps.

  • bin_width_res (int or list of ints) – Width of each residual spectrum bin. If a scalar, the same width is applied to all spectra for all cross spectra. Otherwise, must be a list of up to nspec * nmaps elements, listing bin widths for each of the spectra in res_specs in order, then ordered by map.

  • foreground_fit (bool) – If True, construct bin definitions for foreground components as well.

  • beta_fit (bool) – If True, include delta_beta in the foreground fitting parameters, along with the foreground amplitudes.

  • bin_width_fg (int or list of ints) – Width of each foreground spectrum bin. If a scalar, the same width is applied to all cross spectra. Otherwise, must be a list of up to six elements, listing bin widths for the spectra in the order (TT, EE, BB, TE, EB, TB).

  • lmin_fg (int) – Minimum ell to use for defining foreground bins. If not set, defaults to lmin.

  • lmax_fg (int) – Maximum ell to use for defining foreground bins. If not set, defaults to lmax.

Returns:

bin_def – The bin definition dictionary. Each key contains a Nx2 array that defines the left and right edges for each bin of the corresponding spectrum.

Return type:

dict

get_data_spectra(map_tag=None, transfer=False, do_noise=True)[source]

Return data and noise spectra for the given map tag(s). Data spectra and signal/noise sim spectra must have been precomputed or loaded from disk.

Parameters:
  • map_tag (str) – If None, all map-map cross-spectra are included in the outputs. Otherwise, only the autospectra of the given map are included.

  • transfer (bool) – If True, the data cls are the average of the signal simulations, and noise cls are ignored. If False, the data cls are either cls_data_null (for null tests) or cls_data. See get_masked_data for how these are computed. The input noise is similarly either cls_noise_null or cls_noise.

  • do_noise (bool) – If True, return noise spectra and debiased spectra along with data.

Returns:

  • obs (OrderedDict) – Dictionary of data cross spectra

  • nell (OrderedDict) – Dictionary of noise cross spectra, or None if transfer is True.

  • debias (OrderedDict) – Dictionary of debiased data cross spectra, or None if transfer is True.

get_filename(name, ext='.npz', map_tag=None, iter_index=None, extra_tag=None, data_opts=False, bp_opts=False)[source]

Define a standard output file path to read or write.

Parameters:
  • name (string) – String name of output type. E.g. ‘data_xcorr’ for data cross-correlation spectra. If an output tag is set, the name is appended with ‘_<output_tag>’.

  • ext (string) – File extension. The default (‘.npz’) is used for storing output data dictionaries.

  • map_tag (string) – If supplied, the name is appended with ‘_map_<map_tag>’. Use this argument when storing output data in a loop over input maps.

  • iter_index (int) – If supplied, the name is appended with ‘_iter<iter_index>’.

  • extra_tag (string) – If supplied the extra tag is appended to the name as is.

  • data_opts (bool) – If True, the output filename is constructed by checking the following list of options used in constructing data cross-spectra: ensemble_mean, ensemble_median, sim_index, sim_type, data_type, template_type, reference_type.

  • bp_opts (bool) – If True, also check the following attributes (in addition to those checked if data_opts is True): weighted_bins, return_cls.

Returns:

filename – Output filename as <output_root>/<name><ext>, where <name> is constructed from the above set of options.

Return type:

string

get_files(data_root, data_subset='full/*0', data_root2=None, data_subset2=None)[source]

Find all files for the given data root and subset. Finds all files with data_type 'raw' that match the subset criterion and stores their file tag and frequency information. Subsequent data selection for signal/noise/foreground/template components or other data types should match the set of file tags found here. This method is called at the 'files' checkpoint.

The expected data structure is:

<data_root>
    -> data_<data_type>
        -> full
            -> map_<tag>.fits
            ...
        -> 1of4 (same filenames as full)
        -> 2of4 ('')
        -> 3of4 ('')
        -> 4of4 ('')
    -> signal_<signal_type>
       -> spec_signal_<signal_type>.dat
       -> full
          -> map_<tag>_####.fits
          ...
       -> 1of4 (same filenames as full)
       -> 2of4 (same filenames as full)
       -> 3of4 (same filenames as full)
       -> 4of4 (same filenames as full)
    -> noise_<noise_type> (same filenames as signal_<signal_type>)
    -> masks_<mask_type>
        -> mask_map_<tag>.fits
        ...
    -> foreground_<foreground_type_sim>
        (same filenames as signal_<signal_type>)
    -> templates_<template_type>
       -> template1
          (same filenames as data_<data_type>)
       -> template2
          (same filenames as data_<data_type>)
    -> reference_<reference_type>
       -> reference1
          (same filenames as data_<data_type>)
       -> reference2
          (same filenames as data_<data_type>)
<data_root2> (If provided, same structure as data_root)
    ...

See get_mask_weights, get_masked_data, get_masked_sims for documentation of the various file types.

Parameters:
  • data_root (string) – Top level path containing subdirectories for data, signal sims, noise sims, and masks.

  • data_subset (string) – Subset of maps to use for spectrum estimation. This should be a string that is parseable using glob on the path data_<data_type>/<data_subset>.fits. For example, 'full/*0' will expand to read in the 150 GHz and 90GHz maps. Maps are then sorted in alphabetical order, and identified by their file tag, where each filename is map_<tag>.fits.

  • data_root2 (string) – The root and subset for a second set of data. If either of these is keywords is supplied, then the two data sets are treated as two halves of a null test. In this case, XFaster computes the sum and difference spectra for each map tag in order to estimate a null spectrum.

  • data_subset2 (string) – The root and subset for a second set of data. If either of these is keywords is supplied, then the two data sets are treated as two halves of a null test. In this case, XFaster computes the sum and difference spectra for each map tag in order to estimate a null spectrum.

Returns:

file_settings – A dictionary of file settings used throughout the run. These are stored in full as <output_root>/files_<output_tag>.npz, and a subset are added to the run configuration file <output_root>/config_<output_tag>.txt.

Return type:

dict

get_kernels(window_lmax=None)[source]

Compute kernels using the mask cross-spectra. This follows the polspice azimuthal approximation for the kernel computation.

Parameters:

window_lmax (int) – The window within which the kernel is computed about each ell bin.

Notes

This method is called at the ‘kernels’ checkpoint and loads or saves the following data keys to disk:

kern:

temperature kernel

pkern, mkern:
  • and - polarization kernels

xkern:

temperature/polarization cross term kernel

get_likelihood(qb, inv_fish, map_tag=None, null_first_cmb=False, lmin=33, lmax=250, mcmc=True, alpha_tags='all', beam_tags='all', r_prior=[0, inf], alpha_prior=[0, inf], res_prior=None, beam_prior=[0, 1], num_walkers=50, num_steps=20000, converge_criteria=0.01, reset_backend=None, file_tag=None, r_specs=['ee', 'bb'], template_specs=['ee', 'bb', 'eb'])[source]

Explore the likelihood, optionally with an MCMC sampler.

Parameters:
  • qb (OrderedDict) – Bandpower parameters previously computed by Fisher iteration.

  • inv_fish (array_like) – Inverse Fisher matrix computed with the input qb’s.

  • map_tag (string) – If not None, then the likelihood is sampled using the spectra corresponding to the given map, rather over all possible combinations of map-map cross-spectra. The input qb’s and inv_fish must have been computed with the same option.

  • null_first_cmb (bool) – Keep first CMB bandpowers fixed to input shape (qb=1).

  • lmin (int) – The minimum ell value to be included in the likelihood calculation

  • lmax (int) – The maximum ell value to be included in the likelihood calculation

  • mcmc (bool) – If True, sample the likelihood using an MCMC sampler. Remaining options determine parameter space and sampler configuration.

  • alpha_tags (list of strings) – List of map tags from which foreground template maps should be subtracted. These should be the original map tags, not those generated for chunk sets. If “all”, then all map tags used in the template subtraction are included, as determined by the keys in the template_alpha attribute. If None, then the alpha parameters are not included in the likelihood.

  • beam_tags (list of strings) – List of map tags from which beam error envelopes should be marginalized over. These should be the original map tags, not those generated for chunk sets. If “all”, then all available map tags in the dataset are included. If None, then the beam error parameters are not included in the likelihood.

  • r_prior (2-list or None) – Prior upper and lower bound on tensor to scalar ratio. If None, the fiducial shape spectrum is assumed, and the r parameter space is not varied.

  • alpha_prior (2-list or None) – Prior upper and lower bound on template coefficients. If None, the alpha parameter space is not varied.

  • res_prior (2-list or none) – Prior upper and lower bound on residual qbs. If None, the res parameter space is not varied.

  • beam_prior (2-list or none) – Prior mean and width of gaussian width on beam error (when multiplied by beam error envelope). If None, the beam parameter space is not varied.

  • num_walkers (int) – Number of unique walkers with which to sample the parameter space.

  • num_steps (int) – Maximum number of steps each walker can take in sampling the parameter space.

  • converge_criteria (float) – Convergence criteria for likelihood MCMC chains

  • reset_backend (bool) – If True, clear the backend buffer before sampling. If False, samples are appended to the existing buffer. If not supplied, set to True if the checkpoint has been forced to be rerun.

  • file_tag (string) – If supplied, appended to the likelihood filename.

  • r_specs (list) – Which spectra to use in the r likelihood.

  • template_specs (list) – Which spectra to use for alpha in the likelihood.

get_map(filename, check_nside=True, cache=False, **kwargs)[source]

Load an input map from file or from an internal cache. Maps are checked to make sure they all have a consistent size, and optionally cached to limit disk I/O.

Parameters:
  • filename (string) – Path to file on disk.

  • check_nside (bool) – If True (default), make sure that all maps have the same nside, and that it satisfies lmax <= 4 * nside.

  • cache (bool) – If True, cache the map in memory to avoid rereading from disk. Use this for maps that are used multiple times by the algoritm (e.g. masks).

  • healpy.read_map. (Any remaining arguments are passed to)

Returns:

map – 2D map array containing 1 (T) or 3 (T/Q/U) maps. If the XFaster class was initialized with pol = True, this returns a 2D array of T/Q/U maps from the file. Otherwise a (1, npix) array is returned containing only the T map.

Return type:

array_like

get_mask(filename, cache=True, check_lims=True, **kwargs)[source]

Load an input mask from file or from an internal cache. See XFaster.get_map for details.

Parameters:
  • filename (string) – Path to mask file on disk.

  • cache (bool) – This option defaults to True, since masks are typically used for all data and sims for a given map tag.

  • check_lims (bool) – If True, values in the mask outside of [0,1] are fixed to these limits.

  • XFaster.get_map. (Any remaining arguments are passed to)

Returns:

mask – 2D array containing 1 (T) or 2 (T/P) maps; If the XFaster class was initialized with pol_mask = True, this returns a 2D array containing both T and P masks. Otherwise, a (1, npix) is returned containing only the T map.

Return type:

array_like

get_mask_weights(mask_type='rectangle', apply_gcorr=False, reload_gcorr=False, gcorr_file=None)[source]

Compute cross spectra of the masks for each data map.

Mode counting matrices are also computed and stored for each mask.

Parameters:
  • mask_type (string) – The variant of mask to use, e.g. ‘rectangle’, etc. We assume a mask per file tag in the mask_<mask_type> folder, corresponding to the files in data.

  • apply_gcorr (bool) – If True, a correction factor is applied to the g (mode counting) matrix. The correction factor should have been pre-computed for each map tag using independent scripts in the code package.

  • reload_gcorr (bool) – If True, reload the gcorr file from the masks directory. Useful when iteratively solving for the correction terms.

  • gcorr_file (str) – If not None, path to gcorr file. Otherwise, use file labeled mask_map_<tag>_gcorr.npz in mask directory for signal, or mask_map_<tag>_gcorr_null.npz for nulls.

Notes

This method is called at the ‘masks’ checkpoint, loads or saves a data dictionary with the following keys:

wls:

mask1-x-mask2 mask cross spectra for every mask pair

fsky, w1, w2, w4:

sky fraction and weighted modes per mask product

gmat:

mode-counting matrix, computed from g = fsky * w2 ** 2 / w4

get_masked_data(data_type='raw', template_type=None, template_noise_type=None, template_alpha=None, template_specs=None, reference_type=None, ensemble_mean=False, ensemble_median=False, sim=False, components=['signal', 'noise', 'foreground'], index=None, signal_type_sim=None, r=None, noise_type_sim=None, qb_file=None, foreground_type_sim=None, template_type_sim=None, template_alpha_sim=None, save_sim=False, update_template=False)[source]

Compute cross spectra of the data maps.

If only one dataset is selected, spectra are computed for every combination of pairs of data maps. This results in N * (N + 1) / 2 cross spectra for N maps. A unique mask is used for each input map.

If two datasets are selected for a null test, then sum and difference cross-spectra are computed by summing and differencing the two datasets. A unique mask is used for each map in the first dataset, and the same mask is applied to the corresponding map in the second dataset, so that both halves are masked identically.

If template_type and template_alpha is supplied, the values given are applied to an appropriate template, and the result is subtracted from the data alms with map tags in the dictionary. Map alms are cached to speed up processing, if this method is called repeatedly with different values of template_alpha.

The remaining options handle subtraction of additional biases from the data, or constructing simulated datasets from sim ensembles.

Parameters:
  • data_type (string) – The type of data to use, default: “raw”

  • template_type (string) – Tag for directory (templates_<template_type>) containing templates (e.g. a foreground model) to be scaled by a scalar value per map tag and subtracted from the data. The directory contains one template per map tag.

  • template_noise_type (string) – Tag for directory containing template noise sims to be averaged and scaled similarly to the templates themselves. These averaged sims are used to debias template cross spectra due to correlations in the way the noise ensembles are constructed. Typically, this would be a noise model based on the Planck FFP10 ensemble for each half-mission foreground template. If not supplied, this debiasing step is not performed.

  • template_alpha (dict) – Dictionary of template scaling factors to apply to foreground templates to be subtracted from the data. Keys should match original map tags in the data set.

  • template_specs (list) – Which spectra to use for alpha in the likelihood.

  • reference_type (string) – If supplied, subtract a reobserved reference signal from each data map. The reference signal maps should be two datasets with uncorrelated noise, such as Planck half-mission maps. This option is used for removing expected signal residuals from null tests.

  • ensemble_mean (bool) – If True, use the mean of the signal_type and noise_type ensembles, rather than using maps from the data_type directory or any other sim options. This is useful for testing the behavior of the estimator and mapmaker independently of the data.

  • ensemble_median (bool) – If True, use the median of the signal_type and noise_type ensembles, rather than using maps from the data_type directory or any other sim options. This is useful for testing the behavior of the estimator and mapmaker independently of the data.

  • sim (bool) – If True, construct simulated data maps using the options below.

  • components (list of strings) – List of components to include in the simulated data, of signal, noise, foreground or template.

  • index (int or dict) – If supplied and sim is True, then simulated data maps are constructed from the appropriate index from the sim ensembles signal_type_sim, noise_type_sim and/or foreground_type_sim, rather than using maps from the data_type directory. If an integer, then the same index is used for all ensembles. Otherwise, should be a dictionary keyed by component (signal, tensor, noise, foreground). Additionally, the key default can be used to indicate the index to use for components that are not explicitly enumerated in the dictionary. If not supplied, and sim is True, index 0 is used for all ensembles.

  • signal_type_sim (string) – The variant of signal sims to include in the simulated data maps. This enables having a different signal sim ensemble to use for simulated data than the ensemble from which the signal model component is computed. If this is not supplied and r is not None, then two signal types are searched: ‘signal_r0’ for a scalar component, and ‘signal_r1tens’ for a tensor component. The two maps are linearly combined with a scalar r value to construct a signal map for a simulated dataset.

  • r (float) – If supplied, the simulated signal maps are constructed by combining signal_scalar + r * signal_tensor, where the scalar maps are stored in the directory signal_root_sim, and the tensor maps are stored in the directory tensor_root_sim. A separate sim index is assumed for the signal (scalar) ensemble and the tensor ensemble.

  • noise_type_sim (string) – The variant of noise sims to include in the simulated data maps. This enables having a different noise sim ensemble to use for simulated data than the ensemble from which the noise is computed.

  • qb_file (str) – If supplied and noise is included in components, correct the simulated noise spectra using the noise residuals stored in this file. Typically, this is the output of a separate data run used to determine the appropriate noise correction. See get_noise_residuals for more details.

  • foreground_type_sim (string) – The variant of foreground sims to include in the simulated data maps.

  • template_type_sim (string) – Tag for directory containing foreground templates, to be scaled by a scalar value per map tag and added to the simulated data. The directory contains one template per map tag.

  • template_alpha_sim (dict) – Dictionary of template scaling factors to apply to foreground templates to be added to the simulated the data. Keys should match original map tags in the data set.

  • save_sim (bool) – If True and constructing a simulated dataset using any of the above sim options, write the simulated dataset to disk using an appropriate 'data_xcorr.npz' filename. If False, only non-simulated datasets are written to disk.

  • update_template (bool) – If True, just apply the loaded template with updated template_alpha and template_spec parameters and return.

Notes

This method is called at the ‘data’ checkpoint, loads or saves a data dictionary with the following spectra:

cls_data:

map1-x-map2 cross spectra for every map pair. This contains the sum cross spectra if constructing a null test.

cls_data_clean:

template-subtracted spectra, if template_alpha is supplied.

cls_template:

template cross spectra necessary to rebuild the template-subtracted data when the template_alpha parameter is changed.

cls_data_null:

(map1a-map1b)-x-(map2a-map2b) difference cross spectra for every map pair, if computing a null test

cls_ref, cls_ref_null:

reference cross spectra, if reference_type is set.

get_masked_sims(transfer=False, signal_type='synfast', signal_subset='*', noise_type=None, noise_subset='*', qb_file=None)[source]

Compute average signal and noise spectra for a given ensemble of sim maps. The same procedure that is used for computing data cross spectra is used for each realization in the sim ensemble, and only the average spectra for all realizations are stored.

See get_masked_data for more details on how cross spectra are computed.

Parameters:
  • transfer (bool) – If True, use the ‘sims_transfer’ checkpoint, otherwise use the ‘sims’ checkpoint.

  • signal_type (string) – The variant of signal simulation to use, typically identified by the input spectrum model used to generate it, e.g ‘synfast’. If transfer is True, this signal ensemble is used for the transfer function calculation. The directory may also contain the input Cl spectrum as a camb file, to be used for constructing the kernel-convolved signal model.

  • signal_subset (string) – Subset of map tags to use for spectrum estimation for signal sims. This should be a string that is parseable using glob that is added onto the data_subset path to indicate which sims to use. For example, for all, use '*'. For the first 300 sims, use '0[0-2]*'.

  • noise_type (string) – The variant of noise simulation to use, e.g. ‘gaussian’, ‘stationary’, etc. The directory should contain the same number of simulations for each map tag.

  • noise_subset (string) – Subset of map tags to use for spectrum estimation for noise sims. This should be a string that is parseable using glob that is added onto the data_subset path to indicate which sims to use. For example, for all, use '*'. For the first 300 sims, use '0[0-2]*'.

  • qb_file (string) – Pointer to a bandpowers.npz file in the output directory, used to correct the ensemble mean noise spectrum by the appropriate residual terms. See get_noise_residuals for details.

Notes

This method is called at the ‘sims’ or ‘sims_transfer’ checkpoint, and loads or saves a data dictionary with the following keys:

cls_signal, cls_signal_null:

Mean signal spectra

cls_noise, cls_noise_null:

Mean noise spectra

cls_sim, cls_sim_null:

Mean signal+noise spectra

cls_med, cls_med_null:

Median signal+noise spectra

cls_res, cls_res_null:

NxN, SxN and NxS spectra for computing noise residuals

For null tests, difference spectra of the two null halves are stored in the corresponding *_null keys, and summed spectra are stored in the normal keys. Note that these differ from standard non-null spectra by a factor of 2.

get_masked_template_noise(template_noise_type)[source]

Compute all combinations of template noise cross spectra from the ensemble of template noise sims.

Parameters:

template_noise_type (string) – Tag for directory containing template noise sims to be averaged and scaled similarly to the templates themselves. These averaged sims are used to debias template cross spectra due to correlations in the way the noise ensembles are constructed. Typically, this would be a noise model based on the Planck FFP10 ensemble for each half-mission foreground template.

Returns:

cls_template_noise – Dictionary of template noise spectra averaged over all sims, containing the following keys: [“temp1:temp1”, “temp2:temp2”, “temp1:temp2”]. Each entry has the same shape structure as the cls_data attribute.

Return type:

OrderedDict

Notes

This method is called at the ‘template_noise’ checkpoint, and loads or saves a data dictionary containing the cls_template_noise attribute that is returned.

get_model_spectra(qb, cbl, cls_noise=None, cond_noise=None)[source]

Compute unbinned model spectra from qb amplitudes and a Cbl matrix. Requires pre-loaded bin definitions using get_bin_def or get_transfer.

This method is used internally by fisher_calc.

Parameters:
  • qb (dict of arrays) – Array of bandpowers for every spectrum bin.

  • cbl (dict) – Cbl dict as computed by bin_cl_template.

  • cls_noise (OrderedDict) – If supplied, the noise spectrum is applied to the model spectrum.

  • cond_noise (float) – Conditioning noise amplitude to add to TT, EE and BB autospectra, to improve convergence of the fisher iterations. The noise model is constant cond_noise for EE, BB and 10x that for TT.

Returns:

cls – Model spectra. Keyed by spectrum type, e.g. ‘total_xx’ for the total model spectrom, ‘fg_xx’ for the foreground terms, ‘res_xx’ for the residual (noise) terms, where ‘xx’ is one of the six power spectrum components (tt, ee, bb, te, eb, tb). Each entry in the dictionary is itself a dictionary keyed by map cross, e.g. ‘map1:map1’ for an autospectrum term, ‘map1:map2’ for a cross spectrum, etc, and the map names are the same as those in the map_tags attribute. Each individual spectrum is an array of length lmax + 1.

Return type:

dict of arrays

get_noise_residuals(filename)[source]

Return a dictionary of ell-by-ell noise residual spectra from an output bandpowers file, to be applied to noise Alm’s using healpy.almxfl.

get_signal_shape(filename=None, r=None, component=None, flat=None, filename_fg=None, freq_ref=359.7, beta_ref=1.54, signal_mask=None, transfer=False, save=True)[source]

Load a shape spectrum for input to the Fisher iteration algorithm.

If the spectrum is used as input to get_transfer, it must match the spectrum used to generate the simulations, in order to compute the correct transfer function.

Alternatively, the spectrum can be computed using CAMB for arbitrary values of r, typically used to compute the r likelihood once the bandpowers have been computed.

Finally, the spectrum can be flat in ell^2 Cl. This is typically used as the input shape for computing bandpowers for a null test.

Parameters:
  • filename (string) – Filename for a signal spectrum on disk. If None, and r is None and flat is False, this will search for a spectrum stored in signal_<signal_type>/spec_signal_<signal_type>.dat. Otherwise, if the filename is a relative path and not found, the config directory is searched.

  • r (float) – If supplied and flat is False, a spectrum is computed using CAMB for the given r value. Overrides filename.

  • component ('scalar', 'tensor', 'cmb', 'fg') – If ‘scalar’, and r is not None, return just the r=0 scalar terms in the signal model. If ‘tensor’, return just the tensor component scaled by the input r value. If ‘cmb’ or ‘fg’, return just those signal terms.

  • flat (float) – If given, a spectrum that is flat in ell^2 Cl is returned, with amplitude given by the supplied value. Overrides all other options.

  • filename_fg (string) – Filename for a foreground spectrum on disk. If the filename is a relative path and not found, the config directory is searched. If not supplied, a power law dust model spectrum is assumed.

  • freq_ref (float) – In GHz, reference frequency for dust model. Dust bandpowers output will be at this reference frequency.

  • beta_ref (float) – The spectral index of the dust model. This is a fixed value, with an additive deviation from this value fit for in foreground fitting mode.

  • signal_mask (str array) – Include only these spectra, others set to zero. Options: TT, EE, BB, TE, EB, TB

  • transfer (bool) – If True, this is a transfer function run. If filename is None and r is None and flat is False, will search for a signal spectrum stored in signal_<signal_transfer_type>/spec_signal_<signal_transfer_type>.dat.

  • save (bool) – If True, save signal shape dict to disk.

Returns:

cls – Dictionary keyed by spectrum (cmb_tt, cmb_ee, … , fg_tt, …), each entry containing a vector of length 2 * lmax + 1.

Return type:

OrderedDict

get_transfer(converge_criteria=0.005, iter_max=200, save_iters=False, fix_bb_transfer=False)[source]

Compute the transfer function from signal simulations created using the same spectrum as the input shape.

This raises a ValueError if a negative transfer function amplitude is found.

Parameters:
  • converge_criteria (float) – Maximum fractional change in qb that indicates convergence and stops iteration.

  • iter_max (int) – Maximum number of iterations to perform. if this limit is reached, a warning is issued.

  • save_iters (bool) – If True, the output data from each Fisher iteration are stored in an individual npz file.

  • fix_bb_transfer (bool) – If True, after transfer functions have been calculated, impose the BB xfer is exactly equal to the EE transfer.

Returns:

transfer – Ell-by-ell transfer function for each map

Return type:

OrderedDict

Notes

This method is called at the ‘transfer’ checkpoint, and loads or saves a data dictionary named ‘transfer_all’ with the following entries:

nbins:

number of bins

bin_def:

bin definition dictionary (see get_bin_def)

qb_transfer:

binned transfer function for each map and spectrum component

transfer:

ell-by-ell transfer function for each map and spectrum component

Additionally the final output of fisher_iterate is stored in a dictionary called transfer_map<idx> for each map.

init_log(level='notice', filename=None)[source]

Initialize the logger from the input keyword arguments.

Parameters:
  • level (string, optional, default: "notice") – Verbosity level. Options are “critical”, “error”, “warning”, “notice”, “info”, “debug”, “all”.

  • filename (string, optional) – Logging output filename. Default: None (print to sys.stdout)

kernel_precalc(map_tag=None, transfer=False)[source]

Compute the mixing kernels M_ll’ = K_ll’ * F_l’ * B_l’^2. Called by bin_cl_template to pre-compute kernel terms.

Parameters:
  • map_tag (str) – If supplied, the kernels are computed only for the given map tag (or cross if map_tag is map_tag1:map_tag2). Otherwise, it is computed for all maps and crosses.

  • transfer (bool) – If True, set transfer function to 1 to solve for transfer function qbs.

Returns:

mll – Dictionary of M_ll’ matrices, keyed by component spec and xname.

Return type:

OrderedDict

load_data(name, checkpoint, fields=None, to_attrs=True, shape=None, shape_ref=None, alt_name=None, value_ref=None, optional=None, file_fields=None, **file_opts)[source]

Load xfaster data from an output.npz file on disk.

This method is called throughout the code at various checkpoints. If the data exist on disk, they are loaded and returned. If the data are missing or otherwise incompatible, they are recomputed by the corresponding calling method, and trigger all subsequent data to also be recomputed. Data handling is described in the Notes section for methods that use this functionality.

Parameters:
  • name (string) – The name of the data set. The filename is contructed from this as <output_root>/<name>_<output_tag>.npz. If the file is not found then the data are recomputed.

  • checkpoint (string) – The name of the checkpoint to which this dataset applies. If XFaster is initialized at this checkpoint, or if any of the file checks enabled with the following options fails, all quantities from this point forward are recomputed.

  • fields (list of strings) – List of fields that should be present in the data file. If any are not found, the entire dataset and all subsequent step are recomputed.

  • to_attrs (bool or list of bools or strings) – If True, all items in fields are stored as attributes of the parent object. If A list of booleans, must have the same length as fields; any field for which this list item is True is then stored as an attribute of the object. If any list item is a string, then the corresponding field is stored as an attribute with this new name.

  • shape (tuple of ints) – If set, the field specified by shape_ref is checked to have this shape. If this check fails, then all data are recomputed.

  • shape_ref (string) – The reference field whose shape is checked against shape. If None and shape is set, use the first field in fields.

  • alt_name (string) – Alternative to name argument that will be read if file matching name doesn’t exist.

  • value_ref (dict) – Dictionary of reference values that is checked if simply loading a file from disk instead of recomputing– forces rerun of checkpoints if loaded dictionary differs from value_ref.

  • optional (list of strings) – Fields that, if missing from the data loaded from disk, will not trigger force rerunning of any checkpoints.

  • file_fields (list of strings) – Fields that, if missing from the data loaded from disk, will be searched for in a separate files checkpoint file. This is for backwards compatibility.

  • the (Remaining options are passed to get_filename for constructing)

  • path. (output file)

Returns:

data – If all checks above succeed, the requested data are returned. If any tests fail, None is returned, and all subsequent calls to load_data also return None to trigger recomputing all data that may depend on this dataset. The output dictionary has the additional key ‘output_file’ which is set to the path to the data file on disk.

Return type:

dict

load_map_config(filename)[source]

Load the input map configuration file.

The configuration file should be a file that is readable using ConfigParser. It must contain at least a single section called “frequencies”, with keys for each map tag that may be used by the algorithm. If using the harmonic-domain foreground fitting portions of the algorith, the value of each key should be the observing frequency in GHz that is appropriate for each tag. Otherwise, these frequencies can be any floating point value.

Other optional sections include:

beam: Beam window specifications for each of the tags in “frequencies”. The “beam_product” key should be a path to a .npz file containing a dictionary of beam windows keyed by tag. The “beam_product_error” key should be a path to a similar dictionary containing fraction beam error envelopes, also keyed by tag. See get_beams or get_beam_errors for more details.

fwhm: If using a Gaussian beam model, this section should contain a list of FWHM in arcmin for each such tag in “frequencies”. Keys missing here should be present in the beam product file.

fwhm_err: If using a Gaussian beam model, this section should contain a list of fractional errors on the FWHM for each such tag in “frequencies”. Keys missing here should be present in the beam error product file.

transfer: If present, this section should contain each of the keys in “frequencies”, with the value set to “true” if a transfer function should be computed for the tag, and “false” otherwise (in which case the transfer function will be set to unity for all bins). This option is useful for including, e.g. optimally weighted Planck maps with no transfer function in a joint analysis. If not supplied, it is assumed that a transfer function should be computed for every tag in “frequencies”.

Parameters:

filename (str) – Path to config file

log(message, level=None)[source]

Log a message with the given logging level.

Parameters:
  • message (str) – Log message

  • level (string, default : None) – Logging level. Must be one of “critical”, “error”, “warning”, “notice”, “info”, “debug”, “all”. If not supplied, “all” is assumed.

map2alm(m, pol=None)[source]

Wrapper for healpy.map2alm.

Parameters:
  • m (array_like) – Masked input map for which spherical harmonic alms are computed.

  • pol (bool) – If None, this is set using the value with which the object was initialized.

Returns:

alms – Alms for the input map, computed using the equivalent of healpy.map2alm (m, lmax, pol=self.pol, use_weights=True).

Return type:

array_like

save_config(cfg)[source]

Save a configuration file for the current run on disk. This method is used by xfaster_run to store the config in <output_root>/config_<output_tag>.txt.

Parameters:

cfg (XFasterConfig or RawConfigParser object) – Config object containing all relevant arguments to save to disk.

Returns:

filename – Name of the config file saved to disk.

Return type:

str

save_data(name, from_attrs=[], file_attrs=[], **data)[source]

Save xfaster data to an output .npz file on disk.

Parameters:
  • name (string) – The name of the data set. The filename is contructed from this as <output_root>/<name>_<output_tag>.npz. If the file is not found then the data are recomputed.

  • from_attrs (list of strings) – A list of object attributes which should be stored in the data file.

  • file_attrs (list of strings) – A list of object attributes which contain file paths. The data root is stripped from each path to ensure that only relative paths are written to disk, such that the data root may be changed in the future without triggering rerunning of the entire checkpoint.

  • map_tag (str) – Load the dataset corresponding to this map. See get_filename for documentation.

  • iter_index (int) – Load the dataset corresponding to this iteration index. See get_filename for documentation.

  • bp_opts (bool) – Format output bandpowers file. See get_filename for documentation.

  • extra_tag (str) – Tag to add to file name.

  • dictionary. (Any remaining keyword arguments are added to the output)

Returns:

data – A copy of the data dictionary that was stored to disk. The output dictionary has the additional key ‘output_file’ which is set to the path to the data file on disk.

Return type:

dict

save_state(tag)[source]

Save current object state to disk.

Parameters:

tag (str) – Set the name for the output file to state_<tag>. Otherwise the standard file options are applied to produce the output filename. See get_filename for details.

warn(message)[source]

Log a warning message.

Parameters:

message (str) – Warning log message

class xfaster.xfaster_class.XFasterConfig(defaults=None, default_sec='Uncategorized')[source]

ConfigParser subclass for storing command line options and config.

Class that tracks command-line options for storage to disk.

Parameters:
  • defaults (dict) – Dictionary of overall configuration values. Eg: locals() at beginning of function, or vars(args) from argparse

  • default_sec (string, optional) – The name of the default section in the configuration file.

sort()[source]

Sort the items in each section of the configuration alphabetically.

update(options, section=None)[source]

Update configuration options with a dictionary. Behaves like dict.update() for specified section but also clears options of the same name from the default section.

Parameters:
  • options (dict) – The options to update

  • section (string, optional) – Name of section to update. Default: self.default_sec

write(fp=None, sort=True)[source]

Write an .ini-format representation of the configuration state. Keys are stored alphabetically if sort is True.

Parameters:
  • fp (file object) – If None, write to sys.stdout.

  • sort (bool) – If True, sort items in each section alphabetically.

exception xfaster.xfaster_class.XFasterWarning[source]

Warning generated by the XFaster algorithm.

xfaster_exec

The module containing the master function that performs all the algorithm steps in order, and optionally submits the job to a grid system

class xfaster.xfaster_exec.XFasterJobGroup[source]

Class for parsing xfaster options into a job script, and optionally grouping the jobs together.

add_job(**kwargs)[source]

Add xfaster job to script.

Keyword Arguments:
  • xfaster_run. (Most should correspond to arguments accepted by)

  • to (If job-related arguments are present, they will be passed)

  • set_job_options.

reset()[source]

Initialize to a reset state.

set_job_options(output=None, workdir=None, cput=None, wallt=None, ppn=8, nodes=1, mem=5, env_script=None, omp_threads=None, nice=0, queue=None, job_prefix=None, test=False, pbs=False, dep_afterok=None, exclude=None)[source]

Parse options that control the job script, rather than xfaster. Passed to batch_tools.batch_sub.

Parameters:
  • output (string, optional) – Path for output scheduler files, in a logs subdirectory. If None, use output_root. Overrided by workdir.

  • workdir (string, optional) – If not None, path to output scheduler files. Overrides output.

  • cput (string or float or datetime.timedelta, optional) – Amount of CPU time requested. String values should be in the format HH:MM:SS, e.g. ‘10:00:00’. Numerical values are interpreted as a number of hours.

  • wallt (string or float or datetime.timedelta, optional) – Amount of wall clock time requested. String values should be in the format HH:MM:SS, e.g. ‘10:00:00’. Numerical values are interpreted as a number of hours.

  • ppn (int, optional) – Numper of processes per node

  • nodes (int or string, optional) – Number of nodes to use in job If a string, will be passed as-is to PBS -l node= resource If using SLURM and a string, will overwrite node_list if None

  • mem (float or string, optional) – Amount of memory to request for the job. float values in GB. Or pass a string (eg ‘4gb’) to use directly.

  • env_script (string, optional) – Path to script to source during job script preamble For loading modules, setting environment variables, etc

  • omp_threads (int, optional) – Number of OpenMP threads to use per process

  • nice (int, optional) – Adjust scheduling priority (SLURM only). Range from -5000 (highest priority) to 5000 (lowest priority). Note: actual submitted –nice value is 5000 higher, since negative values require special privilege.

  • queue (string, optional) – The name of the queue to which to submit jobs

  • job_prefix (string, optional) – The name of the job. Default: xfaster

  • test (bool) – If True, only print out the job submission script, don’t submit it.

  • pbs (bool) – If True, use pbs scheduler. Else, use slurm.

  • dep_afterok (string or list of strings) – Dependency. Job ID (or IDs) on which to wait for successful completion, before starting this job

  • exclude (string or list of strings) – List of nodes that will be excluded for job. SLURM-only.

submit(group_by=None, verbose=True, **kwargs)[source]

Submit jobs that have been added.

Parameters:
  • group_by (int, optional) – Group xfaster calls into jobs with this many calls each.

  • verbose (bool, optional) – Print the working directory, and the job ID if submitted successfully.

Returns:

job_ids – The IDs of the submitted jobs

Return type:

list of strings

xfaster.xfaster_exec.xfaster_dump(output_file=None, output_root=None, output_tag=None, checkpoint=None, keys=None, verbose=False)[source]

Print the contents of a set of output archive files.

Parameters:
  • output_file (str) – Path to an output npz archive file.

  • output_root (str) – Path to an xfaster output directory containing npz archive files. Ignored if output_file is set.

  • output_tag (str) – Optional file tag. Ignored if output_file is set.

  • checkpoint (str) – Checkpoint for which to print all matching archive files. Ignored if output_file is set.

  • keys (list of str) – List of keys to print for each matching archive file.

  • verbose (bool) – If True, print the entire contents of the dictionary. Otherwise, truncate each entry to a single line.

xfaster.xfaster_exec.xfaster_main()[source]

Main entry point for command-line interface.

xfaster.xfaster_exec.xfaster_parse(args=None, test=False)[source]

Return a parsed dictionary of arguments for the xfaster execution script.

Parameters:
  • args (list of strings, optional) – If not supplied, read from the command line (sys.argv) by argparse.

  • test (bool, optional) – If True, raise a RuntimeError instead of exiting. Useful for interactive testing.

Returns:

args – Dictionary of parsed options

Return type:

dict

xfaster.xfaster_exec.xfaster_run(config='config_example.ini', output_root=None, output_tag=None, verbose='notice', debug=False, dump_state=False, checkpoint=None, alm_pixel_weights=False, alm_iter=None, add_log=False, data_root='all', data_subset='full/*0', data_root2=None, data_subset2=None, lmin=2, lmax=500, pol=True, pol_mask=True, tbeb=False, bin_width=25, weighted_bins=False, residual_fit=True, bin_width_res=25, res_specs=None, foreground_fit=False, beta_fit=False, bin_width_fg=30, lmin_fg=None, lmax_fg=None, pixwin=True, mask_type='rectangle', window_lmax=None, apply_gcorr=False, reload_gcorr=False, gcorr_file=None, signal_type='synfast', signal_subset='*', signal_transfer_type=None, noise_type='stationary', noise_subset='*', qb_file_sim=None, signal_spec=None, signal_transfer_spec=None, foreground_spec=None, model_r=None, freq_ref=359.7, beta_ref=1.54, data_type='raw', template_type=None, template_noise_type=None, template_alpha_tags=None, template_alpha=None, reference_type=None, ensemble_mean=False, ensemble_median=False, sim_data=False, sim_data_components=['signal', 'noise', 'foreground'], sim_index_signal=None, sim_index_noise=None, sim_index_foreground=None, sim_index_default=0, num_sims=1, signal_type_sim=None, sim_data_r=None, noise_type_sim=None, qb_file_data=None, foreground_type_sim=None, template_type_sim=None, template_alpha_tags_sim=None, template_alpha_sim=None, save_sim_data=False, multi_map=True, bandpower_tag=None, converge_criteria=0.005, cond_noise=1e-05, cond_criteria=5000.0, iter_max=200, save_iters=False, return_cls=False, qb_only=False, fix_bb_transfer=False, null_first_cmb=False, delta_beta_prior=None, like_profiles=False, like_profile_sigma=3.0, like_profile_points=100, likelihood=False, mcmc=True, mcmc_walkers=50, like_converge_criteria=0.01, like_tag=None, like_lmin=33, like_lmax=250, like_r_specs=['EE', 'BB'], like_template_specs=['EE', 'BB', 'EB'], like_alpha_tags='all', alpha_prior=None, r_prior=[-inf, inf], res_prior=None, like_beam_tags='all', beam_prior=None)[source]

Main function for running the XFaster algorithm.

Parameters:
  • config (str) – Configuration file. If path doesn’t exist, assumed to be in xfaster/config/<config>

  • output_root (str) – Directory in which to store output files

  • output_tag (str) – File tag for output files

  • verbose (str) – Logging verbosity level. Can be one of [‘critical’, ‘error’, ‘warning’, ‘notice’, ‘info’, ‘debug’, all’].

  • debug (bool) – Store extra data in output files for debugging.

  • dump_state (bool) – Store current state immediately prior to bandpowers checkpoint. Useful for debugging.

  • checkpoint (str) – If supplied, re-compute all steps of the algorithm from this point forward. Valid checkpoints are [‘files’, ‘masks’, ‘kernels’, ‘sims_transfer’, ‘shape_transfer’, ‘transfer’, ‘sims’, ‘beams’, ‘data’, ‘sim_data’, ‘template_noise’, ‘shape’, ‘bandpowers’, ‘beam_errors’, ‘likelihood’].

  • alm_pixel_weights (bool) – If True, set the use_pixel_weights option to True when computing map Alms using healpy.map2alm. If False, sets the use_weights option to True instead. Note: pixel weights ensure accurate Alm computation, but can only be used for analyses where lmax < 1.5 * nside.

  • alm_iter (int) – If given, set the iter option to the given value when computing map Alms using healpy.map2alm. Using more iterations improves the accuracy of the output Alms. If not set, uses the default number of iterations (3) as defined in healpy. Ignored if alm_pixel_weights is True.

  • add_log (bool) – If True, write log output to a file instead of to STDOUT. The log will be in <output_root>/xfaster-<output_tag>.log. This option is useful for logging to file for jobs that are run directly (rather than submitted).

  • data_root (str) – Root directory where all input files are stored

  • data_subset (str) – The subset of the data maps to include from each data split. Must be a glob-parseable string. Include multiple tags as a comma-delimited sequence enclosed in double quotes.

  • data_root2 (str) – Path for second set of maps for null test. If set, XFaster performs a null test between data_root and data_root2.

  • data_subset2 (str) – The subset of the data maps to include from each data split for the second half of a null split.

  • lmin (int) – Minimum ell at which to start the lowest bin of the output spectra.

  • lmax (int) – Maximum ell for which to compute spectra.

  • pol (bool) – If True, polarization spectra are computed.

  • pol_mask (bool) – If True, a separate mask is applied for Q/U maps.

  • tbeb (bool) – If True, compute TB/EB spectra.

  • bin_width (int or array_like of 6 ints) – Width of each ell bin for each of the six output spectra (TT, EE, BB, TE, EB, TB). EE/BB bins should be the same in order to handle mixing correctly.

  • weighted_bins (bool) – If True, use an lfac-weighted binning operator to construct Cbls. By default, a flat binning operator is used.

  • residual_fit (bool) – If True, include noise residual bins in the estimator.

  • bin_width_res (int or array_like of ints) – Width of each residual spectrum bin. If a scalar, the same width is applied to all spectra for all cross spectra. Otherwise, must be a list of up to nspec * nmaps elements, listing bin widths for each of the spectra in res_specs in order, then ordered by map.

  • res_specs (list of strings) – Spectra to include in noise residual fitting. List values can be any of the cross spectra TT, EE, BB, TE, EB, TB, or EEBB for fitting EE and BB residuals simultaneously. If not supplied, this defaults to EEBB for polarized maps, or TT for unpolarized maps.

  • foreground_fit (bool) – Include foreground residuals in the estimator.

  • beta_fit (bool) – If True, include a fit for delta_beta in the estimator. Otherwise, only fit for foreground amplitudes.

  • bin_width_fg (int or array_like of 6 ints) – Width of each ell bin for each of the six output foreground spectra (TT, EE, BB, TE, EB, TB). EE/BB bins should be the same in order to handle mixing correctly.

  • lmin_fg (int) – Minimum ell to use for defining foreground bins. If not set, defaults to lmin.

  • lmax_fg (int) – Maximum ell to use for defining foreground bins. If not set, defaults to lmax.

  • pixwin (bool) – If True, apply pixel window functions to beam windows.

  • mask_type (str) – The variant of mask to use

  • window_lmax (int) – The size of the window used in computing the mask kernels

  • apply_gcorr (bool) – If True, a correction factor is applied to the g (mode counting) matrix. The correction factor should have been pre-computed for each map tag.

  • reload_gcorr (bool) – If True, reload the gcorr file from the masks directory. Useful when iteratively solving for the correction terms.

  • gcorr_file (str) – If not None, path to gcorr file. Otherwise, use file labeled mask_map_<tag>_gcorr.npy in mask directory for signal, or mask_map_<tag>_gcorr_null.npy for nulls.

  • signal_type (str) – The variant of signal sims to use for the signal component of the covariance model.

  • signal_subset (str) – The subset of the signal sims to include. Must be a glob-parseable string.

  • signal_transfer_type (str) – The variant of signal sims to use for computing the transfer function. If not set, defaults to signal_type.

  • noise_type (str) – The variant of noise sims to use for the noise component of the covariance model.

  • noise_subset (str) – The subset of the noise sims to include. Must be a glob-parseable string.

  • qb_file_sim (str) – If not None, pointer to a bandpowers.npz file in the output directory, to correct the noise ensemble by an appropriate set of residual qb values.

  • signal_spec (str) – The spectrum data file to use for estimating signal component bandpowers. If not supplied, will search for spec_signal_<signal_type>.dat in the signal sim directory.

  • signal_transfer_spec (str) – The spectrum data file used to generate signal sims for computing the signal transfer function. If not supplied, will search for spec_signal_<signal_transfer_type>.dat in the transfer signal sim directory.

  • foreground_spec (string) – The spectrum data file to use for estimating foreground component bandpowers. If not supplied, use a simple power law model for dust foregrounds.

  • model_r (float) – The r value to use to compute a spectrum for estimating bandpowers. Overrides signal_spec.

  • freq_ref (float) – In GHz, reference frequency for dust model. Dust bandpowers output will be at this reference frequency.

  • beta_ref (float) – The spectral index of the dust model. This is a fixed value, with an additive deviation from this value fit for in foreground fitting mode.

  • data_type (str) – The data type to use

  • template_type (str) – Tag for directory (templates_<template_type>) containing templates (e.g. a foreground model) to be scaled by a scalar value per map tag and subtracted from the data. The directory is assumed to contain reference1 and reference2 subdirectories, each containing one template per map tag.

  • template_noise_type (string) – Tag for directory containing template noise sims to be averaged and scaled similarly to the templates themselves. These averaged sims are used to debias template cross spectra due to correlations in the way the noise ensembles are constructed. Typically, this would be a noise model based on the Planck FFP10 ensemble for each half-mission foreground template.

  • template_alpha_tags (list of strings) – List of map tags from which foreground template maps should be subtracted. These should be the original map tags, not those generated for chunk sets.

  • template_alpha (list of floats) – Scalar to be applied to template map for subtraction from each of the data with tags in the list template_alpha_tags.

  • reference_type (string) – Tag for directory containing reobserved reference signals, to be subtracted from each data map. The reference signal maps should be two datasets with uncorrelated noise, such as Planck half-mission maps. This option is used for removing expected signal residuals from null tests.

  • ensemble_mean (bool) – If True, substitute S+N ensemble means for Cls to test for bias in the estimator.

  • ensemble_median (bool) – If True, substitute S+N ensemble median for Cls to test for bias in the estimator.

  • sim_data (bool) – If True, construct simulated data spectra using the options below.

  • sim_data_components (list of strings) – List of components to include in simulated data. May include signal, noise, foreground or template components.

  • sim_index_signal (int) – Sim index to use for the signal component that is included in sim_data_components. If None or < 0, takes the value of sim_index_default.

  • sim_index_noise (int) – Sim index to use for the noise component that is included in sim_data_components. If None or < 0, takes the value of sim_index_default.

  • sim_index_foreground (int) – Sim index to use for the foreground component that is included in sim_data_components. If None or < 0, takes the value of sim_index_default.

  • sim_index_default (int) – Default sim index to use for any component with index < 0 or None in sim_index_<comp>.

  • num_sims (int) – If > 1, repeat the data, bandpowers and likelihood checkpoints this many times, incrementing the value of sim_index_default by 1 each time, starting from the input value. Only used if sim_data is True. All other options remain the same for each iteration.

  • signal_type_sim (str) – The variant of signal sims to use for sim_index data maps. This enables having a different noise sim ensemble to use for sim_index run than the ensemble from which the signal is computed.

  • sim_data_r (float) – If not None, construct the signal component of the simulated data by selecting the appropriate index from an ensemble of scalar and tensor maps, such that the signal component is scalar + r * tensor. This assumes that the tensor simulations are constructed with nt=0, so that the linear relationship holds.

  • noise_type_sim (str) – The variant of noise sims to use for sim_index fake data map. This enables having a different noise sim ensemble to use for sim_index run than the ensemble from which the noise is computed.

  • qb_file_data (str) – If not None, pointer to a bandpowers.npz file in the output directory, to correct the noise component of the simulated data by an appropriate set of residual qb values.

  • foreground_type_sim (str) – Tag for directory (foreground_<foreground_type_sim>) where foreground sims are stored that should be added to the signal and noise sims when running in sim_index mode.

  • template_type_sim (string) – Tag for directory containing foreground templates, to be scaled by a scalar value per map tag and added to the simulated data. The directory contains one template per map tag.

  • template_alpha_tags_sim (list of str) – List of map tags to which foreground template maps should be added, if the template component is included in sim_data_components. These should be the original map tags, not those generated for chunk sets. If None, use the same values as template_alpha_tags.

  • template_alpha_sim (list of floats) – Scalar to be applied to template map for addition to each of the simulated data maps with tags in the list template_alpha_tags. If None, use the same values as template_alpha.

  • save_sim_data (bool) – If True, save data_xcorr file to disk for simulated data.

  • multi_map (bool) – If True, compute all cross-spectra between maps

  • bandpower_tag (str) – Tag to append to bandpowers output file

  • converge_criteria (float) – The maximum fractional change in qb to signal convergence and end iteration

  • cond_noise (float) – The level of regularizing noise to add to EE and BB diagonals.

  • cond_criteria (float) – Threshold on covariance condition number. Above this, regularizing noise will be added to covariance to condition it.

  • iter_max (int) – The maximum number of iterations

  • save_iters (bool) – If True, store the output of each Fisher iteration, in addition to the end result.

  • return_cls (bool) – If True, return C_l spectrum rather than the D_l spectrum

  • qb_only (bool) – If True, do not compute signal window functions or C_l or D_l bandpowers

  • fix_bb_transfer (bool) – If True, after transfer functions have been calculated, impose that the BB transfer function is exactly equal to the EE transfer function.

  • null_first_cmb (bool) – If True, keep first CMB bandpowers fixed to input shape (qb=1).

  • delta_beta_prior (float) – The width of the prior on the additive change from beta_ref. If you don’t want the code to fit for a spectral index different from beta_ref, set this to be None.

  • like_profiles (bool) – If True, compute profile likelihoods for each qb, leaving all others fixed at their maximum likelihood values. Profiles are computed over a range +/–sigma as estimated from the diagonals of the inverse Fisher matrix.

  • like_profile_sigma (float) – Range in units of 1sigma over which to compute profile likelihoods

  • like_profile_points (int) – Number of points to sample along the likelihood profile

  • likelihood (bool) – If True, compute the parameter likelihood

  • mcmc (bool) – If True, sample the parameter likelihood using an MCMC sampler

  • mcmc_walkers (int) – Number of MCMC walkers to use in the likelihood

  • like_converge_criteria (float) – Convergence criteria for likelihood MCMC chains

  • like_tag (str) – Tag to append to likelihood output file

  • like_lmin (int) – The minimum ell value to be included in the likelihood calculation

  • like_lmax (int) – The maximum ell value to be included in the likelihood calculation

  • like_r_specs (list) – Which spectra to use in the r likelihood.

  • like_template_specs (list) – Which spectra to use for alpha in the likelihood.

  • like_alpha_tags (list of strings) – List of map tags from which foreground template maps should be subtracted and fit in the likelihood. If “all”, defaults to template_alpha_tags. If None, alpha fitting in the likelihood is disabled.

  • alpha_prior (list of floats) – Flat prior edges for allowed alpha values in the likelihood. Set to None to not fit for alpha values in the likelihood.

  • r_prior (list of floats) – Flat prior edges for allowed r values in the likelihood.

  • res_prior (list of floats) – Flat prior edges for allowed qb residual values in the likelihood. Set to None to not fit for residual qb values in the likelihood.

  • like_beam_tags (list of strings) – List of map tags from which beam error fields are read in to be fit for in the likelihood. If “all”, then all available map tags in the dataset are included. If None, then beam error fitting in the likelihood is disabled.

  • beam_prior (list of floats) – Gaussian prior mean and number of strandard deviations for beam error. This Gaussian is applied as a prior in fitting for beam error in the likelihood. Set to None to not fit for beam error.

xfaster.xfaster_exec.xfaster_submit(**kwargs)[source]

Submit a single xfaster job. The arguments here should agree exactly with the command line flags for submit mode, with kwargs passed to xfaster_run. Run xfaster --help for help.

spec_tools

More generic functions used by methods in xfaster_class

xfaster.spec_tools.dust_model(ell, pivot=80, amp=34.0, index=-2.28, lfac=True)[source]

Construct a power-law dust power spectrum. Default parameter values are for the Planck LIV best-fit EE dust spectrum.

The functional form is:

model = amp * (ell / pivot) ** index
Parameters:
  • ell (array_like) – Multipoles over which to compute the model

  • pivot (scalar) – Pivot ell at which to refence the amplitude

  • amp (scalar) – Model amplitude

  • index (scalar) – Model spectral index

  • lfac (bool) – If True, multiply Cls by ell*(ell+1)/2/pi

Returns:

model – The dust model computed for each input ell value.

Return type:

array_like

xfaster.spec_tools.get_camb_cl(r, lmax, nt=None, spec='total', lfac=True)[source]

Compute camb spectrum with tensors and lensing.

Parameter values are from arXiv:1807.06209 Table 1 Plik best fit

Parameters:
  • r (float) – Tensor-to-scalar ratio

  • lmax (int) – Maximum ell for which to compute spectra

  • nt (scalar, optional) – Tensor spectral index. If not supplied, assumes slow-roll consistency relation.

  • spec (string, optional) – Spectrum component to return. Can be ‘total’, ‘unlensed_total’, ‘unlensed_scalar’, ‘lensed_scalar’, ‘tensor’, ‘lens_potential’.

  • lfac (bool, optional) – If True, multiply Cls by ell*(ell+1)/2/pi

Returns:

cls – Array of spectra of shape (nspec, lmax + 1). Diagonal ordering (TT, EE, BB, TE).

Return type:

array_like

xfaster.spec_tools.load_camb_cl(filename, lmax=None, pol=None, lfac=True)[source]

Load a CAMB spectrum from a text file. Expects a file with at least two columns, where the first column contains ell values, and the rest are spectrum components (TT, EE, BB, TE). Old-style CAMB files, where the spectra are ordered as (TT, TE, EE, BB), are re-indexed into the new-style ordering.

Parameters:
  • filename (str) – Path to spectrum file on disk.

  • lmax (int) – If supplied, return spectrum up to this ell. Otherwise, all ell’s are included in the output.

  • pol (bool) – If True, include polarization spectra in the output, if present in the file. If False, include only the TT spectrum. Otherwise, include all available spectra.

  • lfac (bool) – If True, multiply Cls by ell*(ell+1)/2/pi

Returns:

cls – Array of spectra of shape (nspec, lmax + 1). Diagonal ordering (TT, EE, BB, TE, EB, TB).

Return type:

array_like

xfaster.spec_tools.scale_dust(nu1, nu2=None, nu0=359.7, beta=1.54, delta=False)[source]

Get the factor by which you must multiply the cross spectrum from maps of frequencies nu1 and nu2 to match the dust power at nu0 given spectral index beta. If nu2 is not given, compute the map-domain factor for nu1 alone. Optionally include terms to construct the derivative with respect to beta if delta=True.

Parameters:
  • nu1 (float) – Frequency of map0 in GHz.

  • nu2 (float) – Frequency of map1 in GHz. If not supplied, return the map-domain multiplicative factor for nu1 only.

  • nu0 (float) – Reference frequency from which to compute relative scaling in GHz.

  • beta (float) – Dust spectral index.

  • delta (bool) – If True, also return the multiplicative beta scaling and its first derivative.

Returns:

  • freq_scale (float) – The relative scaling factor for the dust map or spectrum– multiply by this number to get the dust map or spectrum at the reference frequency

  • beta_scale (float) – The multiplicative scaling factor for beta– multiply the frequency-scaled map or spectrum by beta_scale ** delta_beta to obtain the dust map or spectrum at the adjusted spectral index beta + delta_beta.

  • log_beta_scale (float) – The first derivative of the beta scaling– multiply the frequency-scaled map or spectrum by this number to obtain the derivative of the frequency-scaled map or spectrum with respect to beta.

xfaster.spec_tools.wigner3j(l2, m2, l3, m3)[source]

Wigner 3j symbols computed for all valid values of L, as in:

\[\begin{split}\begin{pmatrix} \ell_2 & \ell_3 & L \\ m_2 & m_3 & 0 \\ \end{pmatrix}\end{split}\]
Parameters:
  • l2 (int) – The ell and m values for which to compute the symbols.

  • m2 (int) – The ell and m values for which to compute the symbols.

  • l3 (int) – The ell and m values for which to compute the symbols.

  • m3 (int) – The ell and m values for which to compute the symbols.

Returns:

  • fj (array_like) – Array of size l2 + l3 + 2, indexed by L

  • lmin (int) – The minimum value of L for which fj is non-zero.

  • lmax (int) – The maximum value of L for which fj is non-zero.