API Reference¶
Command-line Entrypoints¶
runXstack¶
A standalone pipeline code for X-ray spectral shifting and stacking.
X-ray spectral stacking is non-trivial compared to optical spectral stacking. The difficulties arise from two facts:
X-ray has much fewer photon counts (Poisson), meaning that spectral counts and uncertainties cannot be scaled simultaneously (as compared to optical);
X-ray has non-diagonal, complex response, meaning that the response needs to be taken into account when stacking.
To tackle these issues, we develop Xstack: a standalone pipeline code for X-ray spectral stacking. The methodology is to first sum all (rest- frame) PI spectra, without any scaling; and then sum the response files (full response, ARFs*RMFs), with appropriate weighting factors to preserve the overall spectral shape.
The key features of Xstack are:
properly account account for individual spectral contribution to the final stack, by assigning data-driven weighting factors for the responses;
preserve Poisson statistics;
support Galactic absorption correction.
Calling Xstack is simple. For this command line version, it is only a single-line task:
runXstack your_filelist.txt --prefix ./results/stacked_
And you will get:
stacked PI spectrum
./results/stacked_pi.fitsstacked background PI spectrum
./results/stacked_bkgpi.fitsstacked response files:
./results/stacked_arf.fits./results/stacked_rmf.fits
and
./results/stacked_fene.fits, which stores the first contributing energy of each individual source.
Or more sophisticatedly:
runXstack your_filelist.txt --prefix ./results/stacked_ --rsp_weight_method SHP --rsp_project_gamma 2.0 --flux_energy_lo 1.0 --flux_energy_hi 2.3 --nthreads 20 --ene_trc 0.2 --same_rmf AllSourcesUseSameRMF.rmf
If you want to do bootstrap, that is also easy:
runXstack your_filelist.txt --prefix ./results/stacked_ --rsp_weight_method SHP --rsp_project_gamma 2.0 --flux_energy_lo 1.0 --flux_energy_hi 2.3 --nthreads 20 --ene_trc 0.2 --same_rmf AllSourcesUseSameRMF.rmf --resample_method bootstrap --num_bootstrap 100
usage: runXstack [-h] [--prefix PREFIX]
[--rsp_weight_method RSP_WEIGHT_METHOD]
[--rsp_project_gamma RSP_PROJECT_GAMMA]
[--flux_energy_lo FLUX_ENERGY_LO]
[--flux_energy_hi FLUX_ENERGY_HI] [--nthreads NTHREADS]
[--num_bkg_groups NUM_BKG_GROUPS] [--ene_trc ENE_TRC]
[--extended] [--same_rmf SAME_RMF] [--same_target]
[--do_cache] [--bootstrap] [--num_bootstrap NUM_BOOTSTRAP]
[--bootstrap_portion BOOTSTRAP_PORTION]
filelist
Positional Arguments¶
- filelist
text file containing the file names
Named Arguments¶
- --prefix
prefix for output stacked PI, BKGPI, ARF, and RMF files; defaults to
'./results/stacked_'Default:
'./results/stacked_'- --rsp_weight_method
method to calculate RSP weighting factor for each source;
SHP: assuming all sources have same spectral shape (only this mode would require flux_energy_lo and flux_energy_hi),FLX: assuming all sources have same shape and energy flux (weigh by exposure time),LMN: assuming all sources have same shape and luminosity (weigh by exposure/dist^2); defaults toSHPDefault:
'SHP'- --rsp_project_gamma
prior photon index value for projecting RSP matrix onto the output energy channel. This is used in the
SHPmethod, to calculate the weight of each response; defaults to 2.0 (typical for AGN).Default:
2.0- --flux_energy_lo
lower end of the energy range in keV for computing flux, used only when
rsp_weight_method=SHP; defaults to1.0Default:
1.0- --flux_energy_hi
upper end of the energy range in keV for computing flux; used only when
rsp_weight_method=SHP; defaults to2.3Default:
2.3- --nthreads
number of cpus used for RMF shifting
Default:
10- --num_bkg_groups
number of background groups
Default:
10- --ene_trc
energy below which the ARF is manually truncated (e.g., 0.2 keV for eROSITA)
Default:
0.0- --extended
whether or not this is an extended source
Default:
False- --same_rmf
specify the name of common rmf, if all sources are to use the same rmf
- --same_target
stack multiple exposures of the same target (no rest-frame shifting, no Galactic NH correction, direct src/bkg PI summation)
Default:
False- --do_cache
save and load individual rest-frame files
Default:
False- --bootstrap
activate bootstrap mode
Default:
False- --num_bootstrap
number of bootstrap experiments
Default:
10- --bootstrap_portion
portion of sources to resample in each bootstrap experiment
Default:
1.0
Shi-Jiang Chen, Johannes Buchner and Teng Liu (C) 2025 <JohnnyCsj666@gmail.com>
clear_rf_files¶
This is a helper script to clear rest-frame files for a given catalog. The rest-frame files are generated by setting do_cache in Xstack_scripts.Xstack_autoscript.main.
usage: clear_rf_files [-h] [--remove] filelist
Positional Arguments¶
- filelist
text file containing the file names
Named Arguments¶
- --remove
actually deleting individual rest-frame files
Default:
False
Shi-Jiang Chen, Johannes Buchner and Teng Liu (C) 2025 <JohnnyCsj666@gmail.com>
Main Runner¶
Main wrapper module for all spectral shifting+stacking procedures¶
- Authors:
Shi-Jiang Chen (MPE, USTC) Johannes Buchner (MPE) Teng Liu (USTC)
- Email:
- class Xstack.Xstack.XstackRunner(pifile_lst, arffile_lst, rmffile_lst, z_lst, bkgpifile_lst=None, nh_lst=None, srcid_lst=None, rspwt_method='SHP', rspproj_gamma=2.0, int_rng=(1.0, 2.3), sample_rmf=None, sample_arf=None, nh_file=None, Nbkggrp=10, ene_trc=None, extended=False, nthreads=1, bootstrap=False, num_bootstrap=10, bootstrap_portion=1.0, prefix='./results/stacked_', run_cmd=None, same_target=False, do_cache=False)¶
Bases:
objectMain module for X-ray spectral shifting & stacking.
- Parameters:
pifile_lst (list or numpy.ndarray) – The input PI spectrum file list.
arffile_lst (list or numpy.ndarray) – The input ARF file list.
rmffile_lst (list or numpy.ndarray) – The input RMF file list.
z_lst (list or numpy.ndarray) – The redshift list.
bkgpifile_lst (list or numpy.ndarray, optional) – The input background PI spectrum list. Defaults to
None.nh_lst (list or numpy.ndarray, optional) – The Galactic absorption column density list in units of 1 cm^{-2}. Defaults to
None.srcid_lst (list or numpy.ndarray, optional) – The source ID list. Defaults to
None.rspwt_method (str, optional) –
Method for calculating ARFSCAL. Defaults to
SHP. Available methods are:SHP: assuming all sources have same spectral shape, recommendedFLX: assuming all sources have same spectral shape and flux (\(\mathrm{erg\ cm^{-2}\ s^{-1}}\) for point sources; \(\mathrm{erg\ cm^{-2}\ s^{-1}\ deg^{-2}}\) for extended sources).LMN: assuming all sources have same spectral shape and luminosity (\(\mathrm{erg}\ \mathrm{s}^{-1}\) for point sources while \(\mathrm{erg}\ \mathrm{s}^{-1}\ \mathrm{deg}^{-2}\) for extended sources)
rspproj_gamma (float, optional) – The prior photon index value for projecting RSP matrix onto the output energy channel. This is used in the
SHPmethod, to calculate the weight of each response. Defaults to2.0(typical for AGN).int_rng (tuple of (float,float), optional) – The energy (keV) range for computing flux under
SHPmode. Defaults to(1.0, 2.3).sample_rmf (str, optional) – Name of sample RMF. Defaults to
None.sample_arf (str, optional) – Name of sample ARF. Defaults to
None.nh_file (str, optional) –
Galactic absorption profile (absorption factor vs. energy) at 1e20 cm^{-2}. If specified, galactic absorption correction will be applied on the ARF before rest-frame shifting.
Should be in .txt format.
Should also contain the following columns in the first extension:
nhene_ce,nhene_wd,factor.factorshould indicate the absorption factor whennh=1e20.An easy way to obtain the
nh_file: iplottbabs*powerlawwithNh=1e20andPhoIndex=0.0,Norm=1in Xspec.
Nbkggrp (int, optional) – Number of groups with similar background-to-source scaling ratio. Defaults to
10.ene_trc (float, optional) – Truncate energy below which manually set ARF and PI counts to zero. For eROSITA,
ene_trcis typically set as 0.2 keV. Defaults toNone.extended (bool, optional) – Whether or not are the sources to be stacked extended sources. The calculation of response weights would be affected. Defaults to
False, i.e., they are point sources.nthreads (int, optional) – Number of CPUs used in shifting RSP.
bootstrap (bool, optional) – Whether or not to do bootstrap. Defaults to
False.num_bootstrap (int, optional) – Number of bootstrap realizations when
bootstrapmode is activated. Defaults to10.bootstrap_portion (float, optional) – Fraction of sources to be sampled in each bootstrap realization. Defaults to
1.prefix (str, optional) – Prefix for output stacked PI, BKGPI, ARF, and RMF files. Defaults to
./results/stacked_.do_cache (bool, optional) – Whether or not to store and read the intermediate rest-frame files for each individual source. If true, source-wise rest-frame files will be saved under the same path as
pifile. This will save a lot of running time, at the cost of additional disk space (and burden on I/O). Defaults toFalse.run_cmd (str, optional) – Command string recorded into FITS
HISTORYcards for provenance. Defaults toNone.same_target (bool, optional) – Stack multiple exposures of the same target in observed frame. In this mode, rest-frame shifting, Galactic NH correction, and source-to-background scaling are disabled. Source/background PI are summed directly, and full response is weighted using
FLXmode. Defaults toFalse.
Examples
data = XstackRunner( pifile_lst = your_pifile_lst, arffile_lst = your_arffile_lst, rmffile_lst = your_rmffile_lst, z_lst = your_z_lst, bkgpifile_lst = your_bkgpifile_lst, prefix = './results/stacked_', ) data.run() # this will produce the stacked PI, bkgPI, ARF, RMF in one go
- run()¶
Run the main procedure of X-ray spectral shifting and stacking. Do
run_standard_mode()ifsame_target==False, and dorun_same_target_mode()ifsame_target==True.- Returns:
pi_stk_rlz (numpy.ndarray) – Stacked PI spectra from all realizations, with shape
(num_bootstrap, Nchan). Note that ifbootstrap==False, the shape is simply(1, Nchan).pierr_stk_rlz (numpy.ndarray) – Stacked PI err spectra from all realizations, with shape
(num_bootstrap, Nchan). Note that ifbootstrap==False, the shape is simply(1, Nchan).bkgpi_stk_rlz (numpy.ndarray) – Stacked bkg PI spectra from all realizations, with shape
(num_bootstrap, Nchan). Note that ifbootstrap==False, the shape is simply(1, Nchan).bkgpierr_stk_rlz (numpy.ndarray) – Stacked bkg PI err spectra from all realizations, with shape
(num_bootstrap, Nchan). Note that ifbootstrap==False, the shape is simply(1, Nchan).specresp_stk_rlz (numpy.ndarray) – Stacked ARFs from all realizations, with shape
(num_bootstrap, Niene). Note that ifbootstrap==False, the shape is simply(1, Niene).prob_stk_rlz (numpy.ndarray) – Stacked RMF matrices from all realizations, with shape
(num_bootstrap, Niene, Nene). Note that ifbootstrap==False, the shape is simply(1, Niene, Nene).
- run_standard_mode()¶
STANDARD MODE:
In this mode, we will do rest-frame shifting for each source, and then do stacking. Background spectra are scaled to corresponding source region size. Complete treatment of full response. This is the default mode, and is recommended when you are stacking different sources at different redshifts.
- run_same_target_mode()¶
SAME-TARGET MODE:
Stack multiple exposures of the same target in observed frame. In this mode, rest-frame shifting, Galactic NH correction, and source-to-background scaling are disabled. Source/background PI are summed directly, and full response is weighted using
FLXmode.
Utility Modules¶
Module for shifting and stacking PI spectra¶
- Authors:
Shi-Jiang Chen (MPE, USTC) Johannes Buchner (MPE) Teng Liu (USTC)
- Email:
- Xstack.utils.pi.read_pi(pi_fname)¶
Read PI spectrum file.
- Parameters:
pi_fname (str) – Observed-frame pi file to be shifted, in standard OGIP format.
- Returns:
pi_chan (list) – PI channel.
pi_coun (list) – Photon counts in each channel.
z (float) – Redshift if exists.
- Xstack.utils.pi.shift_pi(pi_fname, z, ene_lo=None, ene_hi=None, ene_ce=None, ene_wd=None, rmf_fname=None, ene_trc=None)¶
Shift a single PI to rest-frame.
- Parameters:
pi_fname (str) – Observed-frame pi file to be shifted, in standard OGIP format.
z (float) – Redshift.
ene_lo (numpy.ndarray, optional) – Lower edge of channel energy bin (keV).
ene_hi (numpy.ndarray, optional) – Upper edge of channel energy bin (keV).
ene_ce (numpy.ndarray, optional) – (
ene_lo+ene_hi) / 2ene_wd (numpy.ndarray, optional) – (
ene_hi-ene_lo)rmf_fname (str, optional) – RMF file defining channel-energy conversion, in standard OGIP format. This is optional, unless
ene_loene_hiare not specified.ene_trc (float, optional) – Truncate energy (keV) below which manually set ARF and PI counts to zero. For eROSITA,
ene_trcis typically 0.2 keV.
- Returns:
rest_chan (list) – Rest-frame channel.
rest_coun (list) – Photon counts in each rest-frame channel.
pi_chan (list) – Observed-frame channel.
pi_coun (list) – Photon counts in each observed-frame channel.
- Xstack.utils.pi.calc_pi_error(pi_stk)¶
Round PI counts to integer (so that Poisson applies), and calculate PI uncertainty using Poisson formula.
- Parameters:
pi_stk (numpy.ndarray) – Stacked PI array.
- Returns:
pi_stk (numpy.ndarray) – Stacked PI array. Rounded to nearest integer. E.g., 0.4 -> 0, 0.6 -> 1
pierr_stk (numpy.ndarray) – Stacked PI error array.
- Xstack.utils.pi.calc_bkgpi_error(bkgpi_lst, bkgscal_lst, Nbkggrp=10)¶
Calculate stacked bkg PI counts and uncertainties.
- Parameters:
bkgpi_lst (numpy.ndarray or list) – List of bkg PI spectra.
bkgscal_lst (numpy.ndarray or list) – List of bkg PI scaling factors.
Nbkggrp (int, optional) – Number of background groups with similar
bkgscalto be created. Defaults to10.
- Returns:
pi_stk (numpy.ndarray) – Stacked PI array.
pierr_stk (numpy.ndarray) – Stacked PI error array.
- Xstack.utils.pi.get_bkgscal(src_fname, bkg_fname=None)¶
Get background-to-source scaling ratio, which is calculated as:
(1)¶\[\mathrm{Scaling\ Factor} = \frac{\mathrm{AREASCAL}_{\mathrm{src}}}{\mathrm{AREASCAL}_{\mathrm{bkg}}} \times \frac{\mathrm{BACKSCAL}_{\mathrm{src}}}{\mathrm{BACKSCAL}_{\mathrm{bkg}}} \times \frac{\mathrm{EXPOSURE}_{\mathrm{src}}}{\mathrm{EXPOSURE}_{\mathrm{bkg}}}\]- Parameters:
src_fname (str) – Source PI spectrum name.
bkg_fname (str, optional) – Background PI spectrum name. If not specified, will look for it from the header of
src_fname.
- Returns:
bkgscal – Background-to-source scaling ratio.
- Return type:
float
Notes
Equation
bkgscalapplies to both point sources and extended sources.For eROSITA,
EXPOSUREis the total exposure time during which at least one pixel of the extraction aperture is in the FoV. Since the FoV is scanning over the region during the exposure,EXPOSUREis not the real averaged exposure time per pixel in the region. The real averaged exposure time per pixel, after correcting for such “region-covering incompleteness” issue, is actually:\[T_\mathrm{ave} \equiv \frac{\mathrm{BACKSCL}}{\mathrm{REGAREA}} \times \mathrm{EXPOSURE}\]where
REGAREA/BACKSCALis the region-covering-incompleteness- correcting factor.Note that this is different from Eq. 10 of X. Zhang+2024: the bkg spectra should not only be scaled by
REGAREA, but additionally by the averaged exposure per pixel, which effectively results in an exactly same scaling formula as for the point sources (BACKSCAL*EXPOSURE*AREASCAL).
- Xstack.utils.pi.get_expo(src_fname)¶
Get source exposure time.
- Parameters:
src_fname (str) – Source PI spectrum name.
- Returns:
src_expo – Source exposure time.
- Return type:
float
- Xstack.utils.pi.get_rega(src_fname)¶
Get source geometric area, from the non-standard keyword
REGAREA. For non-eROSITA instrument, return 1.- Parameters:
src_fname (str) – Source PI spectrum name.
- Returns:
src_rega – Source region area (\(\mathrm{deg}^2\)).
- Return type:
float
- Xstack.utils.pi.get_areascal_backscal_corrscal(src_fname)¶
Read scaling keywords from
src_fnameheader.- Parameters:
src_fname (str) – PI spectrum file name.
- Returns:
areascal (float) – Header keyword
AREASCAL. -999 if not found.backscal (float) – Header keyword
BACKSCAL. -999 if not found.corrscal (float) – Header keyword
CORRSCAL. -999 if not found.
- Xstack.utils.pi.make_bkggrpflg(bkgscal_lst, Nbkggrp=10)¶
Group the background spectra into
Ngrpgroups, according to the bkg-to-source scaling ratios.Return an array
bkggrpflg_lstthat tells you which group each background PI spectrum should be assigned to.- Parameters:
bkgscal_lst (list or numpy.ndarray) – List of bkg-to-source scaling ratio (considering both
BACKSCALandEXPOSURE) for each background PI spectrum.Nbkggrp (int, optional) – Number of background groups with similar
bkgscalto be created. Defaults to10.
- Returns:
bkggrpflg_lst (numpy.ndarray) – An array that indicates which group each background PI spectrum should be assigned to (length = len(
bkgscal_lst)).bkgscal_ave_lst (numpy.ndrray) – The average bkg-to-source scaling ratio of each group (length =
Ngrp).
- Xstack.utils.pi.write_pi(chan, pi, pierr=None, pi_fname='stacked_pi.fits', expo=10, rega=1, bkgpi_fname=None, rmf_fname=None, arf_fname=None, spec_type='STACKED', z=None, areascal=1.0, backscal=1.0, corrscal=1.0, run_cmd=None)¶
Write PI spectrum file according to OGIP standards. Assume all spectral files (PI, ARF, RMF) under the same path for
XSPECconvenience.- Parameters:
chan (numpy.ndarray) – Stacked src spectrum channel.
pi (numpy.ndarray) – Stacked src spectrum counts.
pierr (numpy.ndarray, optional) – Stacked src spectrum uncertainty. Defaults to
None(useXSPECPOISSERRby default).pi_fname (str, optional) – Output src spectrum name. Defaults to
stacked_srcpi.fits.expo (int or float, optional) – Stacked exposure. Defaults to
10.rega (int or float, optional) – Stacked region area. Defaults to
1.bkgpi_fname (str, optional) – Stacked bkg PI filename. Defaults to
None.rmf_fname (str, optional) – Stacked RMF filename. Defaults to
None.arf_fname (str, optional) – Stacked ARF filename. Defaults to
None.spec_type (str, optional) –
STACKEDfor stacked spectrum, orRESTFRAMfor individual rest-frame spectrum. Counts are stored as integer inSTACKEDmode, while float inRESTFRAMmode.z (float, optional) – Redshift.
areascal (float, optional) – Header keyword
AREASCAL. Defaults to1.0.backscal (float, optional) – Header keyword
BACKSCAL. Defaults to1.0.corrscal (float, optional) – Header keyword
CORRSCAL. Defaults to1.0.run_cmd (str, optional) – Full command string recorded in
HISTORYfor provenance.
- Return type:
None
- Xstack.utils.pi.write_bkgpi(chan, bkgpi, bkgpierr, bkgpi_fname='stacked_bkgpi.fits', expo=10, rega=1, spec_type='STACKED', z=None, areascal=1.0, backscal=1.0, corrscal=1.0, run_cmd=None)¶
Write bkg PI spectrum file according to OGIP standards. Assume all spectral files (PI, ARF, RMF) under the same path for
XSPECconvenience.- Parameters:
chan (numpy.ndarray) – Stacked bkg spectrum channel.
bkgpi (numpy.ndarray) – Stacked bkg spectrum counts.
bkgpierr (numpy.ndarray) – Stacked bkg spectrum uncertainty.
bkgpi_fname (str, optional) – Output bkg spectrum name. Defaults to
stacked_bkgpi.fits.expo (int or float, optional) – Stacked exposure. Defaults to
10.rega (int or float, optional) – Stacked region area. Defaults to
1.spec_type (str, optional) –
STACKEDfor stacked spectrum, orRESTFRAMfor individual rest-frame spectrum. Counts are stored as float regardless.z (float, optional) – Redshift.
areascal (float, optional) – Header keyword
AREASCAL. Defaults to1.0.backscal (float, optional) – Header keyword
BACKSCAL. Defaults to1.0.corrscal (float, optional) – Header keyword
CORRSCAL. Defaults to1.0.run_cmd (str, optional) – Full command string recorded in
HISTORYfor provenance.
- Return type:
None
- Xstack.utils.pi.make_grpflg(src_fname, grp_fname=None, method='EDGE', rmf_fname='', eelo=None, eehi=None, bkg_fname=None, min_net=0)¶
Add
GROUPINGcolumn to the source PI file.- Parameters:
src_fname (str) – Input source PI file name.
grp_fname (str, optional) – Output grouped PI file name. If not specified, will not create output file.
method (str, optional) –
Grouping method. Available methods:
EDGE: Group by fixed energy bin edges. Edges provided byeeloandeehi.MIN_NET: Group by minimum net counts(src-bkg*bkgscal).Needs to specify the
bkg_fnameandmin_netin each group.
rmf_fname (str, optional) – (for
EDGEmethod) RMF file name. If not specified, the code will automatically search the header ofsrc_fname.eelo (numpy.ndarray, optional) – (for
EDGEmethod) Lower edge of fixed energy bin.eehi (numpy.ndarray, optional) – (for
EDGEmethod) Upper edge of fixed energy bin.bkg_fname (str, optional) – Background file name used in
MIN_NETmode. Defaults toNone. If not specified, will look for it in the header ofsrc_fname.min_net (float or int, optional) – Minimum net counts in each group in
MIN_NETmode. Defaults to0.
- Returns:
grpflg –
GROUPINGcolumn written ingrp_fname.- Return type:
numpy.ndarray
- Xstack.utils.pi.rebin_pi(ene_lo, ene_hi, coun, coun_err, grpflg)¶
Rebin PI file according to
grpflg.- Parameters:
ene_lo (numpy.ndarray) – Lower edge of channel energy bin.
ene_hi (numpy.ndarray) – Upper edge of channel energy bin.
coun (numpy.ndarray) – Photon counts in each channel.
coun_err (numpy.ndarray) – Photon counts error in each channel.
grpflg (numpy.ndarray) – Grouping flag. Must have same length as
ene_loorene_hi.
- Returns:
grpene_lo (numpy.ndarray) – Lower edge of grouped energy bin.
grpene_hi (numpy.ndarray) – Upper edge of grouped energy bin.
grpcoun (numpy.ndarray) – Photon counts in each grouped energy bin.
grpcoun_err (numpy.ndarray) – Photon counts error in each grouped energy bin.
Module for shifting and stacking responses¶
- Authors:
Shi-Jiang Chen (MPE, USTC) Johannes Buchner (MPE) Teng Liu (USTC)
- Email:
- Xstack.utils.rsp.read_rsp(rsp_fname)¶
Read RMF/RSP file.
- Parameters:
rsp_fname (str) – RMF or RSP file name.
- Returns:
prob (numpy.ndarray) – RMF 2D probability matrix, or RSP 2D matrix.
z (float) – Redshift if exists.
- Xstack.utils.rsp.read_rsp_from_arf_rmf(arf_fname, rmf_fname)¶
Read observer-frame full response matrix (RSP = ARF * RMF).
- Parameters:
arf_fname (str) – ARF file name.
rmf_fname (str) – RMF file name.
- Returns:
rspmat – Observer-frame full response matrix.
- Return type:
numpy.ndarray
- Xstack.utils.rsp.shift_rsp(arf_fname, rmf_fname, z, nh_file=None, nh=1e+20, ene_trc=None, ene_lo=None, ene_hi=None)¶
Rest-frame shifting the ARF&RMF. This is literally done by three steps:
Combine GalNH-corrected ARF and RMF into a single RSP matrix (full response);
Shift in the direction of output channel energy. That is to say, shift and broaden the probability profile for each input energy (i.e. when the detector receive a photon with some input energy, the probability that a signal at some output channel energy will be observed; so this is a function of output channel energy) by (1+z);
Shift in the direction of input energy by (1+z), with height (effective area) unchanged.
- Parameters:
arf_fname (str) – The ARF file name.
rmf_fname (str) – The RMF file name.
z (float) – Redshift.
nh_file (str, optional) –
Galactic absorption profile (absorption factor vs. energy). If specified, galactic absorption correction will be applied on the ARF before shifting.
Should be in txt format.
Should also contain the following columns in the first extension:
nhene_ce,nhene_wd,factor.factorshould indicate the absorption factor whennh=1e20.An easy way to obtain the
nh_file: iplottbabs*powerlawwithNh=1e20andPhoIndex=0.0,Norm=1inXSPEC.
nh (float, optional) – The galactic absorption nh of the source (e.g. 3e20). Defaults to
1e20.ene_trc (float, optional) – Truncate energy below which manually set ARF and PI counts to zero. For eROSITA,
ene_trcis typically 0.2 keV. Defaults toNone.
- Returns:
rspmat_sft – Shifted 2D RSP matrix.
- Return type:
numpy.ndarray
- Xstack.utils.rsp.shift_matrix(prob, iene_lo, iene_hi, ene_lo, ene_hi, z)¶
Numba code for Non-parametric RSP/RMF shifting.
- Parameters:
prob (numpy.ndarray) – The RMF 2D probability matrix, or the RSP 2D matrix.
iene_lo (numpy.ndarray) – Lower edge of input model energy (ARF energy) bin.
iene_hi (numpy.ndarray) – Upper edge of input model energy (ARF energy) bin.
ene_lo (numpy.ndarray) – Lower edge of output channel energy bin.
ene_hi (numpy.ndarray) – Upper edge of output channel energy bin.
z (float) – Redshift.
- Returns:
prob_sft – The rest-frame shifted RSP/RMF 2D matrix.
- Return type:
numpy.ndarray
- Xstack.utils.rsp.compute_rspwt(specresp, pi, z, bkgpi, bkgscal, expo, ene_wd, flg, rspwt_method, extended=False, rega=1)¶
Get the weighting factor for a single RSP.
- Parameters:
specresp (numpy.ndarray) – RSP specresp projected on channel energy axis (cm^2 vs. channel energy). This is not simply the ARF curve.
pi (numpy.ndarray) – PI spectrum.
z (float) – Redshift.
bkgpi (numpy.ndarray) – Background PI spectrum.
bkgscal (float) – Background scaling-ratio.
expo (float) – Exposure.
ene_wd (numpy.ndarray) – Output channel energy bin width.
flg (numpy.ndarray) – Output channel energy flag.
method (str) –
Method for calculating ARFSCAL. Available methods are:
SHP: assuming all sources have same spectral shapeFLX: assuming all sources have same spectral shape and fluxFor point sources (
extended==False), flux is in units of erg/cm^2/s.For extended sources (
extended==True), flux is in units of erg/cm^2/s/deg^2.
LMN: assuming all sources have same luminosityFor point sources (
extended==False), luminosity is in units of erg/s.For extended sources (
extended==True), luminosity is in units of erg/s/deg^2.
extended (bool, optional) – Whether or not the source is extended. Defaults to
False, i.e., a point source.rega (int or float, optional) –
REGAREAlist. Used whenextended==True.
- Returns:
rspwt (numpy.ndarray) – The RSP weight for each source.
rspnorm (float) – The RSP weight normalization. This is only useful for
SHPmode.expo_stacked (float) – The final
EXPOSUREto be written in the header of stacked PI and RSP.rega_stacked (float) – The final
REGAREAto be written in the header of stacked PI and RSP.
Notes
The ideal choice of
methodshould beSHP, which starts from the minimum assumption and thus gives the most unbiased results on spectral shape. A caveat ofSHPis that the individual spectrum should have sufficient photon counts (>~10), and the resulting stacked spectrum does not carry a physical flux unit.The second option of
method, in case the individual photon counts is too low, should beFLX. In addition to the minimum assumption used bySHP, it assumes that all sources have similar flux ( erg/cm^2/s for point source or erg/cm^2/s/deg^2 for extended). This should be reasonable for a flux-limited survey, where most sources lie around the detection flux limit, and should thus have similar flux.LMNis similar toFLX, except it assumes all sources to be summed have similar luminosity. Note, luminosity can be calcualted inXSPECas e.g., flux 0.5 2 –> luminosity in 0.5-2 keV band / 1e60
- Xstack.utils.rsp.rescale_rspmat(rspmat, rspwt_lst, expo_lst, rega_lst, rspwt_method, extended=False)¶
Rescale full response matrix (RSP) for different methods.
- Parameters:
rspmat (numpy.ndarray) – Stacked RSP 2D probability matrix.
rspwt_lst (numpy.ndarray) – Response weighting factor for each source (to be rescaled).
expo_lst (numpy.ndarray) – Exposure for each source.
rega_lst (numpy.ndarray) – Region area parameter for each source. TODO: applicable only to eROSITA … update for other inst?
rspwt_method (str) – Response weighting method.
extended (bool, optional) – Extended or not. Defaults to
False.
- Returns:
rspmat (numpy.ndarray) – Rescaled RSP matrix.
rspnorm (float) – To prevent overflow of very large number in the case of
LMNmode, the rescaled RSP matrix has been multiplied by a very small number. Multiply yourrspmatbyrspnormto bring it back to the appropriate number.rspwt_lst (numpy.ndarray) – List of response weighting factors.
expo_stk (float) – Stacked exposure.
rega_stk (float) – Stacked region area.
- Xstack.utils.rsp.correct_arf(specresp, arfene_lo, arfene_hi, factor, nhene_lo, nhene_hi, nh)¶
Multiply the ARF specresp with the galactic absorption profile. The template galactic absorption profile should be at nh=1e20. The source nh value is specified by
nh.- Parameters:
specresp (numpy.ndarray) – The ARF specresp to be corrected.
arfene_lo (numpy.ndarray) – Lower edge of input model energy (ARF energy) bin.
arfene_hi (numpy.ndarray) – Upper edge of input model energy (ARF energy) bin.
factor (numpy.ndarray) – Template galactic absorption profile at
nh=1e20.nhene_lo (numpy.ndarray) – Lower edge of nh model energy bin.
nhene_hi (numpy.ndarray) – Upper edge of nh model energy bin.
nh (float) – Galactic nh of the source.
- Returns:
specresp_cor – The corrected ARF specresp.
- Return type:
numpy.ndarray
- Xstack.utils.rsp.project_rspmat(rspmat, ene_lo, ene_hi, arfene_lo, arfene_hi, proj_axis='CHANNEL', gamma=2.0)¶
Project the 2D RSP matrix onto CHANNEL/MODEL energy axis, to get the effective specresp (cm^2 vs. energy)
- Parameters:
rspmat (numpy.ndarray) – The 2D RSP matrix.
ene_lo (numpy.ndarray) – Lower edge of output channel energy bin.
ene_hi (numpy.ndarray) – Upper edge of output channel energy bin.
arfene_lo (numpy.ndarray) – Lower edge of input model energy (ARF energy) bin.
arfene_hi (numpy.ndarray) – Upper edge of input model energy (ARF energy) bin.
proj_axis (str, optional) –
The projection axis. Available options are:
CHANNEL: project on output channel energy axis. Note that to do this projection, we would nevertheless need to assume a spectral slope, or photo index (specified ingamma). This is to match the convention of unfolded spectrum (in e.g.,XSPEC), where the effective area anchored on channel energy axis is in fact (folded model)/(model).MODEL: project on input model energy axis
Defaults to
CHANNEL.gamma (float, optional) – The spectral slope. Defaults to
2.0. This is only used whenproj_axisisCHANNEL. For AGN sources, a powerlaw with photon index of 2.0 is a good approximation.
- Returns:
rsp1d – The 1D effective area profile.
- Return type:
numpy.ndarray
- Xstack.utils.rsp.get_prob(mat, ebo, f_chan_0=0)¶
Parse the RMF file (input the
MATRIXandEBOUNDSextension) into a 2D probability matrix.- Parameters:
mat (astropy.io.fits.FITS_rec) –
The
MATRIXextension of a standard OGIP RMF file. Must include the following columns:ENERG_LOENERG_HIN_GRPF_CHANN_CHANMATRIX
ebo (astropy.io.fits.FITS_rec) –
The
EBOUNDSextension of a standard OGIP RMF file. Must include the following columns:E_MINE_MAX
f_chan_0 (int, optional) – First channel index. Defaults to
None. If not specified, will be determined from rmf file.
- Returns:
prob – The RMF 2D probability matrix. Index
[i,j], where:irepresents arfene (iene or inpu model energy)jrepresents ene (output channel energy)
- Return type:
numpy.ndarray
- Xstack.utils.rsp.get_prob1d(n_grp, f_chan, n_chan, matrix1d, Nene, f_chan_0=0)¶
Get the 1d probability distribution for output channel energy at a specific input model energy.
- Parameters:
n_grp (int) –
N_GRParray of your specific input model energy, fromMATRIXextension.f_chan (int) –
F_CHANarray of your specific input model energy, fromMATRIXextension.n_chan (int) –
N_CHANarray of your specific input model energy, fromMATRIXextension.matrix1d (numpy.ndarray) –
MATRIXarray of your specific input model energy, fromMATRIXextension.Nene (int) – Length of output channel energy.
f_chan_0 (int, optional) – The index number of the first output channel energy (
0or1). Defaults to0.
- Returns:
prob1d – The 1d probability distribution for output channel energy at a specific input model energy.
- Return type:
numpy.ndarray
- Xstack.utils.rsp.extract_arf_rmf_from_rspmat(rspmat)¶
Extract ARF and RMF from the full response RSP.
- Parameters:
rspmat (numpy.ndarray) – Full response 2D matrix.
- Returns:
specresp (numpy.ndarray) – ARF effective area as a function of input energies (
iene).prob (numpy.ndarray) – 2D probability RMF matrix as a function of input (
iene) and output (ene) energies.
- Xstack.utils.rsp.write_arf(arfene_lo, arfene_hi, specresp, arf_fname='stacked_arf.fits', detchans=1000, expo=10, rega=1, rspwt_method='SHP', rspnorm=1, srcid_lst=None, rspwt_lst=None, pi_totcts_lst=None, bkgpi_totcts_lst=None, flg=None, spec_type='STACKED', z=None, run_cmd=None)¶
Write ARF file according to OGIP standards. Assume all spectral files (PI, ARF, RMF) under the same path for
XSPECconvenience.- Parameters:
arfene_lo (numpy.ndarray) – Lower edge of input model energy (ARF energy) bin, aka
iene_lo.arfene_hi (numpy.ndarray) – Upper edge of input model energy (ARF energy) bin, aka
iene_hi.specresp (numpy.ndarray) – Effective area defined within
arfene_loandarfene_hi.arf_fname (str, optional) – Output ARF name. Defaults to
stacked_arf.fits.detchans (int, optional) – Number of detected channels. This should be the length of PI spectral channels, or equivalently the length of
ene. Defaults to1000.expo (int or float, optional) – Exposure in units of s. Defaults to
10.rega (int or float, optional) – Region area in units of \(\mathrm{deg}^2\). Defaults to
1.rspwt_method (str, optional) – Response weighting method. Defaults to
SHP.rspnorm (int or float, optional) – To prevent overflow of very large number in the case of
LMNmode, the rescaled RSP matrix has been multiplied by a very small number. Multiply yourrspmatbyrspnormto bring it back to the appropriate number. Defaults to1.srcid_lst (numpy.ndarray, optional) – Source id list. Defaults to
None.rspwt_lst (numpy.ndarray, optional) – Response weighting factor list. Defaults to
None.pi_totcts_lst (numpy.ndarray, optional) – PI spectrum total counts. Defaults to
None.bkgpi_totcts_lst (numpy.ndarray, optional) – BKG PI spectrum total counts. Defaults to
None.flg (numpy.ndarray, optional) – Which channels used in
SHPmode in calculating response weighting factors. Defaults toNone.spec_type (str, optional) –
STACKEDif this is the stacked ARF.RESTFRAMif this is single source rest-frame ARF. Defaults toSTACKED.z (float, optional) – Source redshift if this is single source rest-frame ARF.
run_cmd (str, optional) – Full command string recorded in
HISTORYfor provenance.
- Return type:
None
- Xstack.utils.rsp.write_rmf(chan, ene_lo, ene_hi, iene_lo, iene_hi, prob, rmf_fname='./stacked_rmf.fits', expo=10, rega=1, rspwt_method='SHP', srcid_lst=None, rspwt_lst=None, arf_fname='./stacked_arf.fits', spec_type='STACKED', z=None, run_cmd=None)¶
Write RMF file according to OGIP standards. Assume all spectral files (PI, ARF, RMF) under the same path for
XSPECconvenience.- Parameters:
chan (numpy.ndarray) – PI Channel. Must be the same length as
ene.ene_lo (numpy.ndarray) – Lower edge of output channel energy (PI energy) bin.
ene_hi (numpy.ndarray) – Upper edge of output channel energy (PI energy) bin.
prob (numpy.ndarray) – 2D RMF matrix.
rmf_fname (str, optional) – Output RMF name. Defaults to
stacked_rmf.fits.expo (int or float, optional) – Exposure in units of \(\mathrm{s}\). Defaults to
10.rega (int or float, optional) – Region area in units of deg^2. Defaults to
1.rspwt_method (str, optional) – Response weighting method. Defaults to
SHP.srcid_lst (numpy.ndarray, optional) – Source id list. Defaults to
None.rspwt_lst (numpy.ndarray, optional) – Response weighting factor list. Defaults to
None.arf_fname (str, optional) – Associated ARF name. Defaults to
stacked_arf.fits.spec_type (str, optional) –
STACKEDif this is the stacked ARF.RESTFRAMif this is single source rest-frame ARF. Defaults toSTACKED.z (float, optional) – Source redshift if this is single source rest-frame ARF.
run_cmd (str, optional) – Full command string recorded in
HISTORYfor provenance.
- Return type:
None
- Xstack.utils.rsp.get_tlmin_from_header(rmf_fname)¶
Get first channel index from keyword
TLMIN*, according to OGIP standards.- Parameters:
rmf_fname (str) – The RMF file name.
- Returns:
f_chan_0 – First channel index. Will be unity if not found (OGIP default).
- Return type:
int
- Xstack.utils.rsp.rebin_arf(arfene_lo, arfene_hi, specresp, ene_lo, ene_hi, coun, grpflg, prob=None)¶
Anchor the ARF specresp (input model energy) on the output channel energy grid.
- Parameters:
arfene_lo (numpy.ndarray) – Lower edge of input model energy (ARF energy) bin.
arfene_hi (numpy.ndarray) – Upper edge of input model energy (ARF energy) bin.
specresp (numpy.ndarray) – Effective area defined within
arfene_loandarfene_hi.ene_lo (numpy.ndarray) – Lower edge of output channel energy bin.
ene_hi (numpy.ndarray) – Upper edge of output channel energy bin.
coun (numpy.ndarray) – Net photon counts in each channel energy bin.
grpflg (numpy.ndarray) – Channel energy grouping flag, should be passed from
rebin_pi.prob (numpy.ndarray, optional) – The RMF 2D probability matrix. If given, the ARF used for rebinning will be RMF-weighted. Defaults to
None.
- Returns:
grpene_lo (numpy.ndarray) – Lower edge of grouped output channel energy bin.
grpene_hi (numpy.ndarray) – Upper edge of grouped output channel energy bin.
grpspecresp (numpy.ndarray) – Grouped effective area as a function of grouped output channel energy.
- Xstack.utils.rsp.align_arf(ene_lo, ene_hi, arfene_lo, arfene_hi, specresp, prob=None)¶
The ARF energy bin and RMF energy bin (also the PI channel energy bin) does not always match. Align the ARF to get the effective area at each RMF energy bin.
- Parameters:
ene_lo (numpy.ndarray) – Lower edge of output channel energy bin.
ene_hi (numpy.ndarray) – Upper edge of output channel energy bin.
arfene_lo (numpy.ndarray) – Lower edge of input model energy (ARF energy) bin.
arfene_hi (numpy.ndarray) – Upper edge of input model energy (ARF energy) bin.
specresp (numpy.ndarray) – The ARF specresp (cm^2 vs. arf energy).
prob (numpy.ndarray, optional) – RMF 2D matrix (prob.shape=(len(
arfene_lo),len(ene_lo))).
- Returns:
specresp_ali – The aligned ARF specresp.
- Return type:
numpy.ndarray
- Xstack.utils.rsp.concat_rmf(rmf_fname1, rmf_fname2, Es, Ee, Ngrid, out_fname)¶
Concatenate two RMFs into a single large RMF.
- Parameters:
rmf_fname1 (str) – Name of rmf with lower energy.
rmf_fname2 (str) – Name of rmf with higher energy.
Es (float) – Starting energy of the output rmf. Cannot be larger than minimum energy of
rmf_fname1.Ee (float) – Ending energy of the output rmf. Cannot be smaller than maximum energy of
rmf_fname2.Ngrid (int) – Number of grids between
Esandrmf_fname1(also betweenrmf_fname1andrmf_fname2, betweenrmf_fname2andEe).out_fname (str) – Output rmf name.
- Returns:
prob – The output 2D RMF matrix.
- Return type:
numpy.ndarray
- Xstack.utils.rsp.concat_arf(arf_fname1, arf_fname2, Es, Ee, Ngrid, out_fname)¶
Concatenate two ARFs into a single large ARF.
- Parameters:
arf_fname1 (str) – Name of arf with lower energy.
arf_fname2 (str) – Name of arf with higher energy.
Es (float) – Starting energy of the output arf. Cannot be larger than minimum energy of
arf_fname1.Ee (float) – Ending energy of the output arf. Cannot be smaller than maximum energy of
arf_fname2.Ngrid (int) – Number of grids between
Esandarf_fname1(also betweenarf_fname1andarf_fname2, betweenarf_fname2andEe).out_fname (str) – Output ARF name.
- Returns:
specresp – The output ARF specresp.
- Return type:
numpy.ndarray
Module for First energy file (FENE)¶
- Authors:
Shi-Jiang Chen (MPE, USTC) Johannes Buchner (MPE) Teng Liu (USTC)
- Email:
- Xstack.utils.fene.write_fene(srcid_lst, arffene_lst, fene_lst, fene_fname='./stacked_fene.fits', run_cmd=None)¶
Creating a fits storing the first energy of each source’s PI spectrum and ARF specresp.
- Parameters:
srcid_lst (list or numpy.ndarray) – The source ID list.
arffene_lst (list or numpy.ndarray) – The first energy of each sources’s ARF specresp.
fene_lst (list or numpy.ndarray) – The first energy of each source’s PI spectrum.
fene_fname (str) – The output fits name.
run_cmd (str, optional) – Full command string recorded in
HISTORYfor provenance.
- Return type:
None
Miscellaneous tools for Xstack¶
- Authors:
Shi-Jiang Chen (MPE, USTC) Johannes Buchner (MPE) Teng Liu (USTC)
- Email:
- Xstack.utils.misc.get_nh(RA, DEC)¶
Get the Galactic NH from NASA’s HEASARC tool
NH(https://heasarc.gsfc.nasa.gov/Tools/w3nh_help.html). Please ensure the HEASOFT env has been set up. NOTE: this is deprecated. Please usegdpyc.GasMapfrom Github instead.- Parameters:
RA (float)
DEC (float)
- Returns:
nh_val – nh values in units of 1 cm^{-2}
- Return type:
float
- Xstack.utils.misc.pygrppha(src_name, grp_name, grpmin=25)¶
Module for model handling and folding through response¶
- Authors:
Shi-Jiang Chen (MPE, USTC) Johannes Buchner (MPE) Teng Liu (USTC)
- Email:
- Xstack.utils.model.align_model(oarfene_lo, oarfene_hi, omodel, narfene_lo, narfene_hi)¶
Original model (defined on
oarfenegrid) -> New model (defined onnarfenegrid).- Parameters:
oarfene_lo (numpy.ndarray) – Lower edge of original model energy bin.
oarfene_hi (numpy.ndarray) – Upper edge of original model energy bin.
omodel (numpy.ndarray) – Model flux defined on original model energy bin.
narfene_lo (numpy.ndarray) – Lower edge of new model energy bin.
narfene_hi (numpy.ndarray) – Upper edge of new model energy bin.
- Returns:
nmodel – Model flux defined on new model energy bin.
- Return type:
numpy.ndarray
- Xstack.utils.model.fold_model(modelfile, rmffile, arffile, out_name)¶
Fold the input models (\(\mathrm{erg}\ \mathrm{cm}^{-2}\ \mathrm{s}^{-1}\ \mathrm{keV}^{-1}\), input model energy) through response (ARF+RMF) files (\(\mathrm{ct}\ \mathrm{s}^{-1}\ \mathrm{keV}^{-1}\), output channel energy).
Different extensions store different models (models should be defined in
modelfile). Different columns storeE_MIN,E_MAX, and flux of different components in a model.- Parameters:
modelfile (str) – Name of file storing input models to be folded. Different extensions store different models. Different columns store different components.
rmffile (str) – Name of RMF file.
arffile (str) – Name of ARF file.
out_name (str) – Output fits name.
usecpu (int) – Number of CPUs used in folding process.
- Return type:
None
Module for bootstrap & randomization¶
- Authors:
Shi-Jiang Chen (MPE, USTC) Johannes Buchner (MPE) Teng Liu (USTC)
- Email:
- Xstack.utils.random.calc_bootstrap_weights(Nsrc, bootstrap_portion=1.0, rng=None)¶
Generate bootstrap weights (see examples for clarification).
- Parameters:
Nsrc (int) – Number of sources in the list.
bootstrap_portion (float, optional) – The portion of sources to be drawn in each bootstrap realization. Defaults to 1.0, which means drawing the same number of sources as the original list (with replacement).
rng (numpy.random._generator.Generator, optional) – The random number generator to use. If None, a new default generator will be created. Defaults to None.
- Returns:
w – Bootstrap weights list.
- Return type:
numpy.ndarray
Examples
Let’s say we have 20 sources, and we want to do 1 bootstrap realizations.
btw_lst = calc_bootstrap_weights(20)
btw_lstis just how many times each source appear in the current realization:array([0, 0, 0, 2, 2, 0, 1, 1, 2, 1, 3, 0, 1, 2, 1, 0, 0, 0, 2, 2])
And we know that:
0-th source appears 0 time.
1-st source appears 0 time.
2-nd source appears 0 time.
3-rd source appears 2 times.
…
If we want to do 10 bootstrap realizations:
btw_lst_rlz = [ calc_bootstrap_weights(20) for _ in range(10) ] # realizations of bootstrap weight list, (num_bootstrap, Nsrc)
This gives us:
[array([0, 0, 1, 2, 0, 2, 1, 1, 0, 2, 1, 2, 0, 1, 2, 1, 1, 1, 0, 2]), array([0, 1, 1, 1, 2, 4, 1, 0, 0, 2, 2, 3, 2, 0, 0, 0, 0, 1, 0, 0]), array([1, 1, 1, 1, 3, 2, 1, 2, 0, 1, 0, 1, 0, 1, 0, 0, 3, 1, 0, 1]), array([0, 5, 1, 0, 0, 3, 0, 1, 1, 0, 1, 0, 3, 0, 1, 1, 0, 3, 0, 0]), array([1, 1, 2, 0, 1, 1, 2, 0, 0, 1, 0, 2, 1, 2, 1, 2, 0, 3, 0, 0]), array([2, 0, 2, 2, 0, 1, 0, 3, 1, 0, 1, 0, 1, 1, 0, 1, 2, 1, 2, 0]), array([2, 0, 1, 0, 1, 1, 1, 2, 1, 0, 2, 1, 1, 0, 3, 1, 0, 2, 1, 0]), array([0, 1, 2, 0, 1, 0, 0, 1, 2, 0, 1, 1, 1, 2, 1, 1, 1, 1, 2, 2]), array([2, 2, 1, 0, 0, 1, 1, 2, 2, 1, 1, 1, 0, 1, 2, 0, 1, 1, 1, 0]), array([3, 4, 0, 1, 3, 2, 0, 0, 2, 0, 1, 0, 0, 2, 0, 0, 1, 0, 0, 1])]
Module for logging¶
- Authors:
Shi-Jiang Chen (MPE, USTC) Johannes Buchner (MPE) Teng Liu (USTC)
- Email:
- Xstack.utils.logger.get_logger(logname)¶
Define a new customized logger.
- Parameters:
logname (str) – The name of the log file to write to.
- Xstack.utils.logger.utc_now_iso()¶
Get the current UTC time in ISO 8601 format (YYYY-MM-DDTHH:MM:SSZ).
- Xstack.utils.logger.get_ram_gb()¶
Get the current RAM usage of the process in gigabytes (GB).
- Xstack.utils.logger.add_run_cmd_history(header, run_cmd)¶
Add command provenance to FITS
HISTORYcards.
Visualization Helpers¶
Helper functions for visualization.
- Xstack.visual.view.view_rmf(rmf_file, n_grid_i=1000, n_grid=1000, fig=None, ax=None, fig_name=None, cmap='gray_r', log_scale=False, v_min_lbound=1e-06, x_label='Output photon energy (keV)', y_label='Input model energy (keV)')¶
A convenient tool for visualizing 2D RMF matrix. 2D interpolation assumed. You can either call it inside your code to visualize RMF alone side other plots you would like to plot; or you can use this function to produce standalone PNG.
- Parameters:
rmf_file (str) – Name of the RMF file.
n_grid_i (int, optional) – Number of grids for the input model energy (does not have to be the same as the length of
ENERG_LOorENERG_HI). Defaults to1000.n_grid (int, optional) – Number of grids for the output photon energy (does not have to be the same as the length of
E_MINorE_MAX). Defaults to1000.fig (matplotlib.figure.Figure, optional) – The current figure.
ax (matplotlib.axes.Axes, optional) – The current axes.
fig_name (str, optional) – Output figure name. If specified, will create an image.
cmap (str, optional) – cmap. Defaults to
"gray_r".log_scale (bool, optional) – If True, use log-scale for cmap.
v_min_lbound (float, optional) – The lower bound of
v_minfor log-cmap. This means thatLogNorm(vmin=np.max(np.min(prob_new),v_min_lbound), vmax=np.max(prob_new)). Defaults to1e-6.x_label (str, optional) – X label. Defaults to
"Output photon energy (keV)".y_label (str, optional) – Y label. Defaults to
"Input model energy (keV)".
- Returns:
ax (matplotlib.axes.Axes) – The current axes.
im (matplotlib.image.AxesImage) – The image object for the RMF matrix.
cbar (matplotlib.colorbar.Colorbar) – The colorbar object for the RMF matrix.
Examples
from Xstack.visual.view import view_rmf from matplotlib import pyplot as plt fig, ax1 = plt.subplots(1,1,figsize=(4,4)) _, im, cbar = view_rmf("path/to/rmf.fits",fig=fig,ax=ax1,fig_name="rmf.png") #--- optionally modify the colorbar im.set_clim(vmin=1e-5, vmax=1e-1) cbar.update_normal(im) cbar.set_label("Probability",size=14)
- Xstack.visual.view.make_dataarf_plot(src_name, bkg_name=None, arf_name=None, rmf_name=None, grp_name=None, normalize_at=None, outname=None, plot=False, ax=None, **kwargs)¶
Make data/arf plot to visualize the stacked spectral shape.
- Parameters:
src_name (str) – Source spectrum file name.
bkg_name (str, optional) – Background spectrum file name. If not specified, will look for it in the header of
src_name.arf_name (str, optional) – ARF file name. If not specified, will look for it in the header of
src_name.rmf_name (str, optional) – RMF file name. If not specified, will look for it in the header of
src_name.grp_name (str, optional) – Grouping file name. Only uses its
GROUPINGcolumn.normalize_at (int or float, optional) – Output spectrum normalized at some energy (keV). Defaults to None.
outname (str, optional) – Output file name. If not specified, will not create output file name.
plot (bool, optional) – Whether or not to make a plot. Defaults to False.
ax (matplotlib.axes.Axes, optional) – The axes to make the plot. Defaults to None.
**kwargs
- Returns:
grpene_lo (numpy.ndarray) – Lower bounds of grouped energy.
grpene_hi (numpy.ndarray) – Upper bounds of grouped energy.
ratio (numpy.ndarray) – Data/arf ratio, in units of \(\mathrm{ct}\ \mathrm{cm}^{-2} \mathrm{s}^{-1}\ \mathrm{keV}^{-1}\).
ratioerr (numpy.ndarray) – Uncertainty in data/arf ratio.
ax (matplotlib.axes.Axes) – The axes for the plot. None if not specified.
- Xstack.visual.view.valid_energy_range_plot(fene_name, src_name, grp_name, bkg_name, rmf_name, ax=None)¶
Plot:
fraction of sources contributing, and
fraction of net counts (total-background),
as a function of rest-frame energy. These two would facilitate determining a valid energy range for stacked spectrum analysis. Neither the fraction of sources contributing, nor the fraction of net counts can be too low.
- Parameters:
fene_name (str) – The name of the fits file containing first contributing energy. Output of Xstack.
src_name (str) – The source spectrum file name.
grp_name (str) – The grouped source spectrum file to be created.
bkg_name (str) – The background spectrum file name.
rmf_name (str) – The RMF file name.
ax (matplotlib.axes.Axes, optional) – The axes to make the plot. If not specified, will use the current axes. Defaults to None.
- Returns:
ax (matplotlib.axes.Axes) – The axes with the plot.
ax_twinx (matplotlib.axes.Axes) – The twin axes for background fraction plot.
- Xstack.visual.view.gaussian(x, amplitude, mean, stddev)¶
A gaussian function.
- Parameters:
x (float or numpy.ndarray)
amplitude (float)
mean (float)
stddev (float)
- Returns:
pdf – The probability at x.
- Return type:
float or numpy.ndarray
- Xstack.visual.view.get_ene_dsp(ene_ce, prob_lst, fixed_mean=True)¶
Get energy dispersion.
- Parameters:
ene_ce (numpy.ndarray) – Output central channel energy.
prob_lst (numpy.ndarray) – Probability profile for some input model energy (this is a function of output channel energy). Must have same length as
ene_ce.fixed_mean (bool) – If true, the mean energy of the Gaussian will be fixed at the nominal energy (which corresponds to maximal probability).
- Returns:
norm (float) – The Gaussian normalization.
ene_nom (float) – The Gaussian central energy (this is the nominal energy for some input energy).
ene_dsp (float) – The Gaussian width (this is the energy dispersion for some input energy).
- Xstack.visual.view.make_dspmap(mat, ebo, out_name, f_chan_0=None)¶
Make energy dispersion map.
- Parameters:
mat (astropy.io.fits.FITS_rec) – The
MATRIXHDU data from a standard RMF file.ebo (astropy.io.fits.FITS_rec) – The
EBOUNDSHDU data from a standard RMF file.out_name (str) – The output dispersion map name.
f_chan_0 (int, optional) – First channel index. Defaults to None. If not specified, will be determined from rmf file.
- Return type:
None