chandra_aca¶
The chandra_aca package provides functionality related to the Chandra Aspect Camera Assembly and associated aspect functionality.
Plotting¶
-
chandra_aca.plot.
bad_acq_stars
(stars)[source]¶ Return mask of ‘bad’ stars, by evaluating AGASC star parameters.
Parameters: stars – astropy table-compatible set of agasc records of stars. Required fields are [‘CLASS’, ‘ASPQ1’, ‘ASPQ2’, ‘ASPQ3’, ‘VAR’, ‘POS_ERR’] Returns: boolean mask true for ‘bad’ stars
-
chandra_aca.plot.
custom_plt_rcparams
(func)[source]¶ Decorator to make a function use the custom rcParams plt params temporarily. This uses a context manage to ensure original params always get restored.
-
chandra_aca.plot.
plot_compass
(roll)[source]¶ Make a compass plot.
Parameters: roll – Attitude roll for compass plot. Returns: matplotlib figure
-
chandra_aca.plot.
plot_stars
(attitude, catalog=None, stars=None, title=None, starcat_time=None, red_mag_lim=None, quad_bound=True, grid=True, bad_stars=None, plot_keepout=False, ax=None)[source]¶ Plot a catalog, a star field, or both in a matplotlib figure. If supplying a star field, an attitude must also be supplied.
Parameters: - attitude – A Quaternion compatible attitude for the pointing
- catalog – Records describing catalog. Must be astropy table compatible. Required fields are [‘idx’, ‘type’, ‘yang’, ‘zang’, ‘halfw’]
- stars – astropy table compatible set of agasc records of stars Required fields are [‘RA_PMCORR’, ‘DEC_PMCORR’, ‘MAG_ACA’, ‘MAG_ACA_ERR’]. If bad_acq_stars will be called (bad_stars is None), additional required fields [‘CLASS’, ‘ASPQ1’, ‘ASPQ2’, ‘ASPQ3’, ‘VAR’, ‘POS_ERR’] If stars is None, stars will be fetched from the AGASC for the supplied attitude.
- title – string to be used as suptitle for the figure
- starcat_time – DateTime-compatible time. Used in ACASC fetch for proper motion correction. Not used if stars is not None.
- red_mag_lim – faint limit for field star plotting.
- quad_bound – boolean, plot inner quadrant boundaries
- grid – boolean, plot axis grid
- bad_stars – boolean mask on ‘stars’ of those that don’t meet minimum requirements to be selected as acq stars. If None, bad_stars will be set by a call to bad_acq_stars().
- plot_keepout – plot CCD area to be avoided in star selection (default=False)
Returns: matplotlib figure
Star probabilities¶
Functions related to probabilities for star acquisition and guide tracking.
Current default model: grid-floor-2018-11
- The grid-floor-2018-11 model definition and fit values based on:
- https://github.com/sot/aca_stats/blob/master/fit_acq_model-2018-11-binned-poly-binom-floor.ipynb
- See also development models:
- https://github.com/sot/aca_stats/blob/master/fit_acq_model-2018-11-dev
SSAWG initial reviews: 2018-11-07, 2018-11-14 Final review and approval: 2018-11-28
-
chandra_aca.star_probs.
acq_success_prob
(date=None, t_ccd=-10.0, mag=10.0, color=0.6, spoiler=False, halfwidth=120, model=None)[source]¶ Return probability of acquisition success for given date, temperature, star properties and search box size.
Any of the inputs can be scalars or arrays, with the output being the result of the broadcasted dimension of the inputs.
The default probability model is defined by
DEFAULT_MODEL
in this module. Allowed values for themodel
name are ‘sota’, ‘spline’, or a grid model specified by ‘grid-<name>-<date>’ (e.g. ‘grid-floor-2018-11’).Parameters: - date – Date(s) (scalar or np.ndarray, default=NOW)
- t_ccd – CD temperature(s) (degC, scalar or np.ndarray, default=-10C)
- mag – Star magnitude(s) (scalar or np.ndarray, default=10.0)
- color – Star color(s) (scalar or np.ndarray, default=0.6)
- spoiler – Star spoiled (boolean or np.ndarray, default=False)
- halfwidth – Search box halfwidth (arcsec, default=120)
- model – probability model name: ‘sota’ | ‘spline’ | ‘grid-*’
Returns: Acquisition success probability(s)
-
chandra_aca.star_probs.
binom_ppf
(k, n, conf, n_sample=1000)[source]¶ Compute percent point function (inverse of CDF) for binomial, where the percentage is with respect to the “p” (binomial probability) parameter not the “k” parameter. The latter is what one gets with scipy.stats.binom.ppf().
This function internally generates
n_sample
samples of the binomal PMF in order to compute the CDF and finally interpolate to get the PPF.The following example returns the 1-sigma (0.17 - 0.84) confidence interval on the true binomial probability for an experiment with 4 successes in 5 trials.
Example:
>>> binom_ppf(4, 5, [0.17, 0.84]) array([ 0.55463945, 0.87748177])
Parameters: - k – int, number of successes (0 < k <= n)
- n – int, number of trials
- conf – float or array of floats, percent point values
- n_sample – number of PMF samples for interpolation
Returns: percent point function values corresponding to
conf
-
chandra_aca.star_probs.
clip_and_warn
(name, val, val_lo, val_hi, model)[source]¶ Clip
val
to be in the rangeval_lo
toval_hi
and issue a warning if clipping occurs. Thename
andmodel
are just used in the warning.Parameters: - name – Value name
- val – Value
- val_lo – Minimum
- val_hi – Maximum
- model – Model name
Returns: Clipped value
-
chandra_aca.star_probs.
get_box_delta
(halfwidth)[source]¶ Transform from halfwidth (arcsec) to the box_delta value which gets added to failure probability (in probit space).
Parameters: halfwidth – scalar or ndarray of box sizes (halfwidth, arcsec) Returns: box deltas
-
chandra_aca.star_probs.
grid_model_acq_prob
(mag=10.0, t_ccd=-12.0, color=0.6, halfwidth=120, probit=False, model=None)[source]¶ Calculate a grid model probability of acquisition success for a star with specified mag, t_ccd, color, and search box halfwidth.
This does a 3-d linear interpolation on mag, t_ccd, and halfwidth using a pre-computed gridded model that is stored in a FITS file.
Parameters: - mag – ACA magnitude (float or np.ndarray)
- t_ccd – CCD temperature (degC, float or ndarray)
- color – B-V color to check for B-V=1.5 => red star (float or np.ndarray)
- halfwidth – search box halfwidth (arcsec, default=120, float or ndarray)
- probit – if True then return Probit(p_success). Default=False
- model – Model name, e.g. ‘grid-floor-2018-11’
Returns: Acquisition success probability(s)
-
chandra_aca.star_probs.
guide_count
(mags, t_ccd, count_9th=False)[source]¶ Calculate a guide star fractional count/metric using signal-to-noise scaled mag thresholds. This uses a modification of the guide star fractional counts that were suggested at the 7-Mar-2018 SSAWG and agreed upon at the 21-Mar-2018 SSAWG. The implementation here does a piecewise linear interpolation between the reference mag - fractional count points instead of the original “threshold interpolation” (nearest neighbor mag <= reference mag). Approved at 16-Jan-2019 SSAWG. One feature is the slight incline in the guide_count curve from 1.0005 at mag=6.0 to 1.0 at mag=10.0. This does not show up in standard outputs of guide_counts to two decimal places (8 * 0.0005 = 0.004), but helps with minimization. :param mags: magnitude(s) :param t_ccds: CCD temperature(s) :param count_9th: return fractional count of 9th mag or brighter stars :returns: fractional count
-
chandra_aca.star_probs.
mag_for_p_acq
(p_acq, date=None, t_ccd=-10.0, halfwidth=120, model=None)[source]¶ For a given
date
andt_ccd
, find the star magnitude that has an acquisition probability ofp_acq
. Star magnitude is defined/limited to the range 5.0 - 12.0 mag.Parameters: - p_acq – acquisition probability (0 to 1.0)
- date – observation date (any Chandra.Time valid format)
- t_ccd – ACA CCD temperature (deg C)
- halfwidth – search box halfwidth (arcsec, default=120)
- model – probability model (see acq_success_prob() for allowed values, default)
Returns mag: star magnitude
-
chandra_aca.star_probs.
model_acq_success_prob
(mag, warm_frac, color=0, halfwidth=120)[source]¶ Call sota_model_acq_prob() with same params. This is retained purely for back-compatibility but use is deprecated.
-
chandra_aca.star_probs.
prob_n_acq
(star_probs)[source]¶ Given an input array of star acquisition probabilities
star_probs
, return the probabilities of acquiring exactly n_acq stars, where n_acq is evaluated at values 0 to n_stars. This is returned as an array of length n_stars. In addition the cumulative sum, which represents the probability of acquiring n_acq or fewer stars, is returned.Parameters: star_probs – array of star acq probabilities (list or ndarray) Returns n_acq_probs, cum_n_acq_probs: tuple of ndarray, ndarray
-
chandra_aca.star_probs.
set_acq_model_ms_filter
(ms_enabled=False)[source]¶ Choose one of two sets of acquisition model fit parameters based on
ms_enabled
. The default isFalse
:- True: MS filtering enabled (prior to FEB0816 loads), where stars would be rejected if MS flag was set
- False: MS filtering disabled (including and after FEB0816 loads)
The selected fit parameters are global/module-wide.
-
chandra_aca.star_probs.
sota_model_acq_prob
(mag, warm_frac, color=0, halfwidth=120)[source]¶ Calculate raw SOTA model probability of acquisition success for a star with
mag
magnitude and a CCD warm fractionwarm_frac
. This is not typically used directly since it does not account for star properties like spoiled or color=0.7.Uses the empirical relation:
P_fail_probit = offset(mag) + scale(mag) * warm_frac + box_delta(halfwidth) P_acq_fail = Normal_CDF() P_acq_success = 1 - P_acq_fail
This is based on the dark model and acquisition success model presented in the State of the ACA 2013, and subsequently updated to use a Probit transform and separately fit B-V=1.5 stars. It was updated in 2017 to include a fitted dependence on the search box halfwidth. See:
https://github.com/sot/skanb/blob/master/pea-test-set/fit_box_size_acq_prob.ipynb https://github.com/sot/aca_stats/blob/master/fit_acq_prob_model-2017-07-sota.ipynb
Parameters: - mag – ACA magnitude (float or np.ndarray)
- warm_frac – N100 warm fraction (float or np.ndarray)
- color – B-V color to check for B-V=1.5 => red star (float or np.ndarray)
- halfwidth – search box halfwidth (arcsec, default=120)
Returns: Acquisition success probability(s)
-
chandra_aca.star_probs.
spline_model_acq_prob
(mag=10.0, t_ccd=-12.0, color=0.6, halfwidth=120, probit=False)[source]¶ Calculate poly-spline-tccd model (aka ‘spline’ model) probability of acquisition success for a star with specified mag, t_ccd, color, and search box halfwidth.
The model definition and fit values based on: - https://github.com/sot/aca_stats/blob/master/fit_acq_prob_model-2018-04-poly-spline-tccd.ipynb
See also: - Description of the motivation and initial model development.
- Final review and approval.
- https://occweb.cfa.harvard.edu/twiki/bin/view/Aspect/StarWorkingGroupMeeting2018x04x18
Parameters: - mag – ACA magnitude (float or np.ndarray)
- t_ccd – CCD temperature (degC, float or ndarray)
- color – B-V color to check for B-V=1.5 => red star (float or np.ndarray)
- halfwidth – search box halfwidth (arcsec, default=120, float or ndarray)
- probit – if True then return Probit(p_success). Default=False
Returns: Acquisition success probability(s)
-
chandra_aca.star_probs.
t_ccd_warm_limit
(mags, date=None, colors=0, halfwidths=120, min_n_acq=5.0, cold_t_ccd=-16, warm_t_ccd=-5, model=None)[source]¶ Find the warmest CCD temperature which meets the
min_n_acq
acquisition stars criterion. This returns a value betweencold_t_ccd
andwarm_t_ccd
. At the cold end the result may be belowmin_n_acq
, in which case the star catalog may be rejected.The
min_n_acq
argument can be specified in one of two ways:- Scalar float value: expected number of acquired stars must exceed threshold.
- Expected number is the sum of the individual probabilities.
- Tuple (n, prob): computed probability of acquiring
n
or fewer stars - must not exceed
prob
.
- Tuple (n, prob): computed probability of acquiring
Parameters: - mags – list of star ACA mags
- date – observation date (any Chandra.Time valid format)
- colors – list of star B-V colors (optional, default=0.0)
- halfwidths – list of acq box halfwidths(optional, default=120)
- min_n_acq – float or tuple (see above)
- cold_t_ccd – coldest CCD temperature to consider (default=-16 C)
- warm_t_ccd – warmest CCD temperature to consider (default=-5 C)
- model – probability model (see acq_success_prob() for allowed values, default)
Returns: (t_ccd, n_acq | prob_n_or_fewer) tuple with CCD temperature upper limit and: - number of expected ACQ stars at that temperature (scalar min_n_acq) - probability of acquiring
n
or fewer stars (tuple min_n_acq)
-
chandra_aca.star_probs.
t_ccd_warm_limit_for_guide
(mags, min_guide_count=4.0, warm_t_ccd=-5.0, cold_t_ccd=-16.0)[source]¶ Solve for the warmest temperature that still gets the min_guide_count. This returns a value between
cold_t_ccd
andwarm_t_ccd
. At the cold end the result may be belowmin_n_acq
, in which case the star catalog may be rejected.Parameters: - mags – list of star ACA mags
- min_guide_count – float minimum fractional guide count
- warm_t_ccd – warmest CCD temperature to consider (default=-5 C)
- cold_t_ccd – coldest CCD temperature to consider (default=-16 C)
Returns: t_ccd
Transformations¶
The transform modules includes:
- Coordinate conversions between ACA Pixels and ACA Y-angle, Z-angle
- Mag to count conversions
- Science target coordinate to ACA frame conversions
-
chandra_aca.transform.
broadcast_arrays
(*args)[source]¶ Broadcast *args inputs to same shape and return an
is_scalar
flag and the broadcasted version of inputs. This lets intermediate code work on arrays that are guaranteed to be the same shape and at least a 1-d array, but reshape the output at the end.Parameters: args – tuple of scalar / array inputs Returns: [is_scalar, *flat_args]
-
chandra_aca.transform.
broadcast_arrays_flatten
(*args)[source]¶ Broadcast *args inputs to same shape and then return that shape and the flattened view of all the inputs. This lets intermediate code work on all scalars or all arrays that are the same-length 1-d array and then reshape the output at the end (if necessary).
Parameters: args – tuple of scalar / array inputs Returns: [shape, *flat_args]
-
chandra_aca.transform.
calc_aca_from_targ
(targ, y_off, z_off, si_align=None)[source]¶ Calculate the PCAD (ACA) pointing attitude from target attitude and Y,Z offset (per OR-list).
Parameters: - targ – science target from OR/Obscat (Quat-compliant object)
- y_off – Y offset (deg, sign per OR-list convention)
- z_off – Z offset (deg, sign per OR-list convention)
- si_align – SI ALIGN matrix (default=ODB_SI_ALIGN)
Return type: q_aca (Quat)
-
chandra_aca.transform.
calc_targ_from_aca
(aca, y_off, z_off, si_align=None)[source]¶ Calculate the target attitude from ACA (PCAD) pointing attitude and Y,Z offset (per OR-list).
Parameters: - aca – ACA (PCAD) pointing attitude (any Quat-compliant object)
- y_off – Y offset (deg, sign per OR-list convention)
- z_off – Z offset (deg, sign per OR-list convention)
- si_align – SI ALIGN matrix
Return type: q_targ (Quat)
-
chandra_aca.transform.
count_rate_to_mag
(count_rate)[source]¶ Convert count rate in e- / sec to ACA mag
Based on $CALDB/data/chandra/pcad/ccd_char/acapD1998-12-14ccdN0002.fits columns mag0=10.32 and cnt_rate_mag0=5263.0.
Parameters: count_rate – count rate (e-/sec) Returns: magnitude (ACA mag)
-
chandra_aca.transform.
mag_to_count_rate
(mag)[source]¶ Convert ACA mag to count rate in e- / sec
Based on $CALDB/data/chandra/pcad/ccd_char/acapD1998-12-14ccdN0002.fits columns mag0=10.32 and cnt_rate_mag0=5263.0. To convert to raw telemetered counts (ADU), divide e- by the gain of 5.0 e-/ADU.
Parameters: mag – star magnitude in ACA mag Returns: count rate (e-/sec)
-
chandra_aca.transform.
pixels_to_yagzag
(row, col, allow_bad=False, flight=False, t_aca=20, pix_zero_loc='edge')[source]¶ Convert ACA row/column positions to ACA y-angle, z-angle. It is expected that the row and column input arguments have the same length.
The
pix_zero_loc
parameter controls whether the input pixel values are assumed to have integral values at the pixel center or at the pixel lower/left edge. The PEA flight coefficients assume lower/left edge, and that is the default. However, it is generally more convenient when doing centroids and other manipulations to use the center.Parameters: - row – ACA pixel row (single value, list, or 1-d numpy array)
- col – ACA pixel column (single value, list, or 1-d numpy array)
- allow_bad – boolean switch. If true, method will not throw errors if the row/col values are nominally off the ACA CCD.
- flight – Use flight EEPROM coefficients instead of default ground values.
- t_aca – ACA temperature (degC) for use with flight (default=20C)
- pix_zero_loc – row/col coords are integral at ‘edge’ or ‘center’
Return type: (yang, zang) each vector of the same length as row/col
-
chandra_aca.transform.
radec_to_yagzag
(ra, dec, q_att)[source]¶ Given RA, Dec, and pointing quaternion, determine ACA Y-ang, Z-ang. The input
ra
anddec
values can be 1-d arrays in which case the outputyag
andzag
will be corresponding arrays of the same length.This is a wrapper around Ska.quatutil.radec2yagzag but uses arcsec instead of deg for yag, zag.
Parameters: - ra – Right Ascension (degrees)
- dec – Declination (degrees)
- q_att – ACA pointing quaternion
Returns: yag, zag (arcsec)
-
chandra_aca.transform.
snr_mag_for_t_ccd
(t_ccd, ref_mag, ref_t_ccd, scale_4c=None)[source]¶ Given a t_ccd, solve for the magnitude that has the same expected signal to noise as ref_mag / ref_t_ccd.
If scale_4c is None, the value from dark_model.DARK_SCALE_4C is used.
To solve for the magnitude that has the expected signal to noise as ref_mag / ref_t_ccd, we define this equality:
counts(mag) / noise(t_ccd) = counts(ref_mag) / noise(ref_t_ccd)
We then assume (with some handwaving) that:
noise(t_ccd) = noise(ref_t_ccd) * scale
where
scale = scale_4c ** ((t_ccd - ref_t_ccd) / 4.0)
And we use the definition of counts as:
counts(mag) = ACA_CNT_RATE_MAG0 * 10.0 ** ((ACA_MAG0 - mag) / 2.5)
Substituting in, reducing, and solving the original equality for mag gives
ref_mag - (t_ccd - ref_t_ccd) * np.log10(scale_4c) / 1.6
-
chandra_aca.transform.
yagzag_to_pixels
(yang, zang, allow_bad=False, pix_zero_loc='edge')[source]¶ Convert ACA y-angle/z-angle positions to ACA pixel row, column. It is expected that the y-angle/z-angle input arguments have the same length.
The
pix_zero_loc
parameter controls whether the input pixel values are assumed to have integral values at the pixel center or at the pixel lower/left edge. The PEA flight coefficients assume lower/left edge, and that is the default. However, it is generally more convenient when doing centroids and other manipulations to usepix_zero_loc='center'
.Parameters: - yang – ACA y-angle (single value, list, or 1-d numpy array)
- zang – ACA z-angle (single value, list, or 1-d numpy array)
- allow_bad – boolean switch. If true, method will not throw errors if the resulting row/col values are nominally off the ACA CCD.
- pix_zero_loc – row/col coords are integral at ‘edge’ or ‘center’
Return type: (row, col) each vector of the same length as row/col
-
chandra_aca.transform.
yagzag_to_radec
(yag, zag, q_att)[source]¶ Given ACA Y-ang, Z-ang and pointing quaternion determine RA, Dec. The input
yag
andzag
values can be 1-d arrays in which case the outputra
anddec
will be corresponding arrays of the same length.This is a wrapper around Ska.quatutil.yagzag2radec but uses arcsec instead of deg for yag, zag.
Parameters: - yag – ACA Y angle (arcsec)
- zag – ACA Z angle (arcsec)
- q_att – ACA pointing quaternion
Returns: ra, dec (arcsec)
Centroid residuals¶
-
class
chandra_aca.centroid_resid.
CentroidResiduals
(start, stop)[source]¶ Class to calculate star centroid residuals.
This class is designed to set up and perform the residual calculations on any desired combination of source centroids and source attitudes. For the common use cases, centroids, attitudes, and commanded star positions are retrieved automatically from archived sources.
Based on analysis, time offsets are applied to centroid times by default. See fit notebooks in:
http://nbviewer.jupyter.org/url/cxc.harvard.edu/mta/ASPECT/ipynb/centroid_time_offsets/OR.ipynb
and
http://nbviewer.jupyter.org/url/cxc.harvard.edu/mta/ASPECT/ipynb/centroid_time_offsets/ER.ipynb
Users should see the class method
for_slot
for a convenient way to get centroid residuals on anobsid
for an ACAslot
(aka image number).Example usage:
>>> import numpy as np >>> from chandra_aca.centroid_resid import CentroidResiduals >>> cr = CentroidResiduals.for_slot(obsid=20001, slot=5) >>> np.max(np.abs(cr.dyags)) 0.87602233734844503 >>> np.max(np.abs(cr.dzags)) 1.2035827855862777 >>> cr.atts[0] array([-0.07933254, 0.87065874, -0.47833673, 0.08278696]) >>> cr.agasc_id 649201816
This example calculates the residuals on slot 5 of obsid 20001 using the ground aspect solution and ground centroids. Here is another example that does the same thing without using the
for_slot
convenience.Example usage:
>>> import numpy as np >>> from chandra_aca.centroid_resid import CentroidResiduals >>> cr = CentroidResiduals(start='2017:169:18:54:50.138', stop='2017:170:05:13:58.190') >>> cr.set_atts('ground') >>> cr.set_centroids('ground', slot=5) >>> cr.set_star(agasc_id=649201816) >>> cr.calc_residuals() >>> np.max(np.abs(cr.dyags)) 0.87602233734844503
Parameters: - start – start time of interval for residuals (DateTime compatible)
- stop – stop time of interval for residuals (DateTime compatible)
-
calc_residuals
()[source]¶ Calculate residuals based on attitude and ra/dec of star. Note that the sampling and times of yags may be different from zags so these should be done independently.
Residuals are available in self.dyags and self.dzags. Predicted values from attitude and star position in self.pred_yags and self.pred_zags
-
set_atts
(source)[source]¶ Get attitude solution quaternions from
source
.One could also just set atts and att_times attributes directly.
-
set_centroids
(source, slot, alg=8, apply_dt=True)[source]¶ Assign centroids from
source
andslot
to the objects centroid attributes (yag, zag, yag_times, zag_times)For the supported sources (ground, obc) the centroids are fetched from the mica L1 archive or telemetry.
yag, zag, yag_times an zag_times can also be set directly without use of this method.
Parameters: - source – ‘ground’ | ‘obc’
- slot – ACA slot
- alg – for ground processing, use centroids from this algorithm.
- apply_dt – apply centroid time offsets via ‘set_offsets’
-
set_offsets
()[source]¶ Apply time offsets to centroids based on type and source of centroid, obsid (suggesting 8x8 or 6x6 data), telemetry source (‘maude’ or ‘cxc’) and aspect solution source. These time offsets were fit. See fit notebooks at:
http://nbviewer.jupyter.org/url/cxc.harvard.edu/mta/ASPECT/ipynb/centroid_time_offsets/OR.ipynb
and
http://nbviewer.jupyter.org/url/cxc.harvard.edu/mta/ASPECT/ipynb/centroid_time_offsets/ER.ipynb
ACA alignment drift¶
Function(s) related to ACA alignment.
In particular compute the dynamical pointing offset required to achieve a desired zero-offset target aimpoint.
A key element of this module is the fitting analysis here: https://github.com/sot/aimpoint_mon/blob/master/fit_aimpoint_drift-2018-11.ipynb
-
class
chandra_aca.drift.
AcaDriftModel
(scale, offset, trend, jumps, year0)[source]¶ Define a drift model for aspect solution SIM DY/DZ values as a function of time and ACA CCD temperature. This expresses the model which is defined and fitted in the fit_aimpoint_drift notebook in this repo.
-
calc
(times, t_ccd)[source]¶ Calculate the drift model for aspect solution SIM DY/DZ values for input
times
andt_ccd
. The two arrays are broadcasted to match.The returned drifts are in arcsec and provide the expected aspect solution SIM DY or DZ values in mm. This can be converted to a drift in arcsec via the scale factor 20.493 arcsec/mm.
Parameters: - times – array of times (CXC secs)
- t_ccd – CCD temperatures (degC)
Returns: array of ASOL SIM DY/DZ (mm)
-
-
chandra_aca.drift.
get_aca_offsets
(detector, chip_id, chipx, chipy, time, t_ccd)[source]¶ Compute the dynamical ACA offset values for the provided inputs.
The
time
andt_ccd
inputs can be either scalars or arrays.Parameters: - detector – one of ACIS-I, ACIS-S, HRC-I, HRC-S
- chipx – zero-offset aimpoint CHIPX
- chipy – zero-offset aimpoint CHIPY
- chip_id – zero-offset aimpoint CHIP ID
- time – time(s) of observation (any Chandra.Time compatible format)
- t_ccd – ACA CCD temperature(s) (degC)
Returns: aca_offset_y, aca_offset_z (arcsec)
-
chandra_aca.drift.
get_default_zero_offset_table
()[source]¶ Get official SOT MP zero offset aimpoint table:
Returns: zero offset aimpoint table as astropy.Table
-
chandra_aca.drift.
get_target_aimpoint
(date, cycle, detector, too=False, zero_offset_table=None)[source]¶ Given date, proposal cycle, and detector, return aimpoint chipx, chipy, chip_id
Parameters: - date – observation date
- cycle – proposal cycle of observation
- detector – target detector
- too – boolean. If target is TOO use current cycle not proposal cycle.
- zero_offset_able – table (astropy or numpy) of zero offset aimpoint table defaults to official SOT MP version if not supplied.
Returns: tuple of chipx, chipy, chip_id
ACA images¶
The aca_image module contains classes and utilities related to ACA readout images. This includes the ACAImage class for manipulating images in ACA coordinates, a first moment centroiding routine, and a library for generating synthetic ACA images that have a high-fidelity point spread function (PSF).
-
class
chandra_aca.aca_image.
ACAImage
[source]¶ ACAImage is an ndarray subclass that supports functionality for the Chandra ACA. Most importantly it allows image indexing and slicing in absolute “aca” coordinates, where the image lower left coordinate is specified by object
row0
andcol0
attributes.It also provides a
meta
dict that can be used to store additional useful information. Any keys which are all upper-case will be exposed as object attributes, e.g.img.BGDAVG
<=>img.meta['BGDAVG']
. Therow0
attribute is a proxy forimg.meta['IMGROW0']
, and likewise forcol0
.When initializing an
ACAImage
, additional*args
and**kwargs
are used to try initializing vianp.array(*args, **kwargs)
. If this fails thennp.zeros(*args, **kwargs)
is tried. In this way one can either initialize from array data or create a new array of zeros.Examples:
>>> import numpy as np >>> from chandra_aca.aca_image import ACAImage >>> dat = np.random.uniform(size=(1024, 1024)) >>> a = ACAImage(dat, row0=-512, col0=-512) >>> a = ACAImage([[1,2], [3,4]], meta={'BGDAVG': 5.2}) >>> a = ACAImage(shape=(1024, 1024), row0=-512, col0=-512)
Parameters: - row0 – row coordinate of lower left image pixel (int, default=0)
- col0 – col coordinate of lower left image pixel (int, default=0)
- meta – dict of object attributes
- *args – additional args passed to np.array() or np.zeros()
- **kwargs – additional kwargs passed to np.array() or np.zeros()
-
aca
¶ Return a light copy (same data) of self but with the _aca_coords attribute switched on so that indexing is absolute.
-
centroid_fm
(bgd=None, pix_zero_loc='center', norm_clip=None)[source]¶ First moment centroid of
self
using 6x6 mousebitten image for input 6x6 or 8x8 images.Note that the returned
norm
is the sum of the background-subtracted 6x6 mousebitten image, not the entire image.Parameters: - bgd – background to subtract, scalar or NxN ndarray (float)
- pix_zero_loc – row/col coords are integral at ‘edge’ or ‘center’
- norm_clip – clip image norm at this min value (default is None and implies Exception for non-positive norm)
Returns: row, col, norm float
-
flicker_init
(flicker_mean_time=10000, flicker_scale=1.0, seed=None)[source]¶ Initialize instance variables to allow for flickering pixel updates.
The
flicker_scale
can be interpreted as follows: if the pixel was going to flicker by a multiplicative factor of (1 + x), now make it flicker by (1 + x * flicker_scale). This applies for flickers that increase the amplitude. For flickers that make the value smaller, then it would be 1 / (1 + x) => 1 / (1 + x * flicker_scale).The flicker_cdf file here was created using: /proj/sot/ska/www/ASPECT/ipynb/chandra_aca/flickering-pixel-model.ipynb
Examples and performance details at: /proj/sot/ska/www/ASPECT/ipynb/chandra_aca/flickering-implementation.ipynb
The model was reviewed and approved at SS&AWG on 2019-05-22.
Parameters: - flicker_mean_time – mean flickering time (sec, default=10000)
- flicker_scale – multiplicative factor beyond model default for flickering amplitude (default=1.0)
- seed – random seed for reproducibility (default=None => no seed)
-
flicker_update
(dt, use_numba=True)[source]¶ Propagate the image forward by
dt
seconds and update any pixels that have flickered during that interval.This has the option to use one of two implementations. The default is to use the numba-based version which is about 6 times faster. The vectorized version is left in for reference.
Parameters: - dt – time (secs) to propagate image
- use_numba – use the numba version of updating (default=True)
-
chandra_aca.aca_image.
centroid_fm
(img, bgd=None, pix_zero_loc='center', norm_clip=None)[source]¶ First moment centroid of
img
.Return FM centroid in coords where lower left pixel of image has value (0.0, 0.0) at the center (for pix_zero_loc=’center’) or the lower-left edge (for pix_zero_loc=’edge’).
Parameters: - img – NxN ndarray
- bgd – background to subtract, float of NXN ndarray
- pix_zero_loc – row/col coords are integral at ‘edge’ or ‘center’
- norm_clip – clip image norm at this min value (default is None and implies Exception for non-positive norm)
Returns: row, col, norm float
-
class
chandra_aca.aca_image.
AcaPsfLibrary
(filename=None)[source]¶ Access the ACA PSF library, whch is a library of 8x8 images providing the integrated (pixelated) ACA PSF over a grid of subpixel locations.
Example:
>>> from chandra_aca.aca_image import AcaPsfLibrary >>> apl = AcaPsfLibrary() # Reads in PSF library data file >>> img = apl.get_psf_image(row=-10.456, col=250.123, norm=100000) >>> img <ACAImage row0=-14 col0=247 array([[ 39, 54, 56, 52, 37, 33, 30, 21], [ 79, 144, 260, 252, 156, 86, 67, 36], [ 162, 544, 2474, 5269, 2012, 443, 163, 57], [ 255, 1420, 10083, 12688, 11273, 1627, 302, 78], [ 186, 1423, 8926, 8480, 12292, 2142, 231, 64], [ 80, 344, 1384, 6509, 4187, 665, 111, 43], [ 40, 78, 241, 828, 616, 188, 65, 29], [ 24, 39, 86, 157, 139, 69, 48, 32]])>
Parameters: filename – file name of ACA PSF library (default=built-in file) Returns: AcaPsfLibrary object -
get_psf_image
(row, col, norm=1.0, pix_zero_loc='center', interpolation='bilinear', aca_image=True)[source]¶ Get interpolated ACA PSF image that corresponds to pixel location
row
,col
.Parameters: - row – (float) row value of PSF centroid
- col – (float) col value of PSF centroid
- norm – (float) summed intensity of PSF image
- pix_zero_loc – row/col coords are integral at ‘edge’ or ‘center’
- interpolation – ‘nearest’ | ‘bilinear’ (default)
- aca_image – return ACAImage if True, else return ndarray
Returns: ACAImage if (aca_image is True) else (ndarray image, row0, col0)
-
ACAImage class¶
ACAImage is an ndarray subclass that supports functionality for the Chandra
ACA. Most importantly it allows image indexing and slicing in absolute
“aca” coordinates, where the image lower left coordinate is specified
by object row0
and col0
attributes.
It also provides a meta
dict that can be used to store additional useful
information. Any keys which are all upper-case will be exposed as object
attributes, e.g. img.BGDAVG
<=> img.meta['BGDAVG']
. The row0
attribute is a proxy for img.meta['IMGROW0']
, and likewise for col0
.
Creation¶
When initializing an ACAImage
, additional *args
and **kwargs
are
used to try initializing via np.array(*args, **kwargs)
. If this fails
then np.zeros(*args, **kwargs)
is tried. In this way one can either
initialize from array data or create a new array of zeros.
One can easily create a new ACAImage
as shown below. Note that it must
always be 2-dimensional. The initial row0
and col0
values default
to zero.
>>> from chandra_aca.aca_image import ACAImage
>>> im4 = np.arange(16).reshape(4, 4)
>>> a = ACAImage(im4, row0=10, col0=20)
>>> a
<ACAImage row0=10 col0=20
array([[ 0, 1, 2, 3],
[ 4, 5, 6, 7],
[ 8, 9, 10, 11],
[12, 13, 14, 15]])>
One could also initialize by providing a meta
dict:
>>> a = ACAImage(im4, meta={'IMGROW0': 10, 'IMGCOL0': 20, 'BGDAVG': 5.2})
>>> a
<ACAImage row0=10 col0=20
array([[ 0, 1, 2, 3],
[ 4, 5, 6, 7],
[ 8, 9, 10, 11],
[12, 13, 14, 15]])>
Note
The printed representation of an ACAImage
is always shown as the
rounded integer version of the values, but the full float value
is stored internally. If you print()
the image then the floating
point values are shown.
Image access and setting¶
You can access array elements as usual, and in fact do any normal numpy array operations:
>>> a[3, 3]
15
The special nature of ACAImage
comes by doing array access via the aca
attribute.
In this case all index values are in absolute coordinates, which might be negative. In
this case we can access the pixel at ACA coordinates row=13, col=23, which is equal to
a[3, 3]
for the given row0
and col0
offset:
>>> a.aca[13, 23]
15
Creating a new array by slicing adjusts the row0
and col0
values
like you would expect:
>>> a2 = a.aca[12:, 22:]
>>> a2
<ACAImage row0=12 col0=22
array([[500, 11],
[ 14, 15]])>
You can set values in absolute coordinates:
>>> a.aca[11:13, 21:23] = 500
>>> a
<ACAImage row0=10 col0=20
array([[ 0, 1, 2, 3],
[ 4, 500, 500, 7],
[ 8, 500, 500, 11],
[ 12, 13, 14, 15]])>
Now let’s make an image that represents the full ACA CCD and set a
sub-image from our 4x4 image a
. This uses the absolute location
of a
to define a slice into b
:
>>> b = ACAImage(shape=(1024,1024), row0=-512, col0=-512)
>>> b[a] = a
>>> b.aca[8:16, 18:26]
<ACAImage row0=8 col0=18
array([[ 0, 0, 0, 0, 0, 0, 0, 0],
[ 0, 0, 0, 0, 0, 0, 0, 0],
[ 0, 0, 0, 1, 2, 3, 0, 0],
[ 0, 0, 4, 500, 500, 7, 0, 0],
[ 0, 0, 8, 500, 500, 11, 0, 0],
[ 0, 0, 12, 13, 14, 15, 0, 0],
[ 0, 0, 0, 0, 0, 0, 0, 0],
[ 0, 0, 0, 0, 0, 0, 0, 0]])>
You can also do things like adding 100 to every pixel in b
within the area of a
:
>>> b[a] += 100
>>> b.aca[8:16, 18:26]
<ACAImage row0=8 col0=18
array([[ 0, 0, 0, 0, 0, 0, 0, 0],
[ 0, 0, 0, 0, 0, 0, 0, 0],
[ 0, 0, 100, 101, 102, 103, 0, 0],
[ 0, 0, 104, 600, 600, 107, 0, 0],
[ 0, 0, 108, 600, 600, 111, 0, 0],
[ 0, 0, 112, 113, 114, 115, 0, 0],
[ 0, 0, 0, 0, 0, 0, 0, 0],
[ 0, 0, 0, 0, 0, 0, 0, 0]])>
Image arithmetic operations¶
In addition to doing image arithmetic operations using explicit slices
as shown previously, one can also use normal arithmetic operators like
+
(for addition) or +=
for in-place addition.
When the right-side operand is included via its ``.aca`` attribute, then the operation is done in ACA coordinates.
This means that the operation is only done on overlapping pixels. This is shown in the examples below. The supported operations for this are:
- Addition (
+
and+=
) - Subtraction (
-
and-=
) - Multiplication (
*
and*=
) - Division (
/
and/=
) - True division (
/
and/=
in Py3+ and with __future__ division) - Floor division (
//
and//=
) - Modulus (
%
and%=
) - Power (
**
and**=
)
Initialize images (different shape and offset)
>>> a = ACAImage(shape=(6, 6), row0=10, col0=20) + 1
>>> a
<ACAImage row0=10 col0=20
array([[1, 1, 1, 1, 1, 1],
[1, 1, 1, 1, 1, 1],
[1, 1, 1, 1, 1, 1],
[1, 1, 1, 1, 1, 1],
[1, 1, 1, 1, 1, 1],
[1, 1, 1, 1, 1, 1]])>
>>> b = ACAImage(np.arange(1, 17).reshape(4, 4), row0=8, col0=18) * 10
>>> b
<ACAImage row0=8 col0=18
array([[ 10, 20, 30, 40],
[ 50, 60, 70, 80],
[ 90, 100, 110, 120],
[130, 140, 150, 160]])>
Add images (output has shape and row0/col0 of left side input)
>>> a + b.aca
<ACAImage row0=10 col0=20
array([[111, 121, 1, 1, 1, 1],
[151, 161, 1, 1, 1, 1],
[ 1, 1, 1, 1, 1, 1],
[ 1, 1, 1, 1, 1, 1],
[ 1, 1, 1, 1, 1, 1],
[ 1, 1, 1, 1, 1, 1]])>
>>> b + a.aca
<ACAImage row0=8 col0=18
array([[ 10, 20, 30, 40],
[ 50, 60, 70, 80],
[ 90, 100, 111, 121],
[130, 140, 151, 161]])>
>>> b += a.aca
>>> b
<ACAImage row0=8 col0=18
array([[ 10, 20, 30, 40],
[ 50, 60, 70, 80],
[ 90, 100, 111, 121],
[130, 140, 151, 161]])>
Make ``b`` image be fully contained in ``a``
>>> b.row0 = 11
>>> b.col0 = 21
>>> a += b.aca
>>> a
<ACAImage row0=10 col0=20
array([[ 1, 1, 1, 1, 1, 1],
[ 1, 11, 21, 31, 41, 1],
[ 1, 51, 61, 71, 81, 1],
[ 1, 91, 101, 112, 122, 1],
[ 1, 131, 141, 152, 162, 1],
[ 1, 1, 1, 1, 1, 1]])>
Normal image addition fails if shape is mismatched
>>> a + b
Traceback (most recent call last):
File "<ipython-input-19-f96fb8f649b6>", line 1, in <module>
a + b
File "chandra_aca/aca_image.py", line 68, in _operator
out = op(self, other) # returns self for inplace ops
ValueError: operands could not be broadcast together with shapes (6,6) (4,4)
Meta-data¶
Finally, the ACAImage
object can store arbitrary metadata in the
meta
dict attribute. However, in order to make this convenient and
distinct from native numpy attributes, the meta
attributes should
have UPPER CASE names. In this case they can be directly accessed
as object attributes instead of going through the meta
dict:
>>> a.IMGROW0
10
>>> a.meta
{'IMGCOL0': 20, 'IMGROW0': 10}
>>> a.NEWATTR = 'hello'
>>> a.meta
{'IMGCOL0': 20, 'NEWATTR': 'hello', 'IMGROW0': 10}
>>> a.NEWATTR
'hello'
>>> a.meta['fail'] = 1
>>> a.fail
Traceback (most recent call last):
AttributeError: 'ACAImage' object has no attribute 'fail'
Maude Decom¶
Classes and functions to help fetching ACA telemetry data using Maude. These include the following global variables
- PIXEL_MAP: dict of np.array, with values mapping integer pixel indices to pixel string ID
- PIXEL_MAP_INV: dict of dict, with values mapping pixel string ID to integer pixel indices.
- PIXEL_MASK: dict of np.array. Values are boolean masks that apply to images of different sizes
- ACA_MSID_LIST: dictionary of commonly-used ACA telemetry MSIDs.
- ACA_SLOT_MSID_LIST: dictionary of ACA image telemetry MSIDs.
PIXEL_MAP contains maps between pixel indices and pixel IDm depending on the image size:
- Size 4X41:
-----------------------------------------
| -- | -- | -- | -- | -- | -- | -- | -- |
-----------------------------------------
| -- | -- | -- | -- | -- | -- | -- | -- |
-----------------------------------------
| -- | -- | D1 | H1 | L1 | P1 | -- | -- |
-----------------------------------------
| -- | -- | C1 | G1 | K1 | O1 | -- | -- |
-----------------------------------------
| -- | -- | B1 | F1 | J1 | N1 | -- | -- |
-----------------------------------------
| -- | -- | A1 | E1 | I1 | M1 | -- | -- |
-----------------------------------------
| -- | -- | -- | -- | -- | -- | -- | -- |
-----------------------------------------
| -- | -- | -- | -- | -- | -- | -- | -- |
-----------------------------------------
- Size 6X61 or 6X62:
-----------------------------------------
| -- | -- | -- | -- | -- | -- | -- | -- |
-----------------------------------------
| -- | -- | E2 | F2 | G2 | H2 | -- | -- |
-----------------------------------------
| -- | D2 | D1 | H1 | L1 | P1 | I2 | -- |
-----------------------------------------
| -- | C2 | C1 | G1 | K1 | O1 | J2 | -- |
-----------------------------------------
| -- | B2 | B1 | F1 | J1 | N1 | K2 | -- |
-----------------------------------------
| -- | A2 | A1 | E1 | I1 | M1 | L2 | -- |
-----------------------------------------
| -- | -- | P2 | O2 | N2 | M2 | -- | -- |
-----------------------------------------
| -- | -- | -- | -- | -- | -- | -- | -- |
-----------------------------------------
- Size 8X81, 8X82, 8X83 or 8X84:
-----------------------------------------
| H1 | P1 | H2 | P2 | H3 | P3 | H4 | P4 |
-----------------------------------------
| G1 | O1 | G2 | O2 | G3 | O3 | G4 | O4 |
-----------------------------------------
| F1 | N1 | F2 | N2 | F3 | N3 | F4 | N4 |
-----------------------------------------
| E1 | M1 | E2 | M2 | E3 | M3 | E4 | M4 |
-----------------------------------------
| D1 | L1 | D2 | L2 | D3 | L3 | D4 | L4 |
-----------------------------------------
| C1 | K1 | C2 | K2 | C3 | K3 | C4 | K4 |
-----------------------------------------
| B1 | J1 | B2 | J2 | B3 | J3 | B4 | J4 |
-----------------------------------------
| A1 | I1 | A2 | I2 | A3 | I3 | A4 | I4 |
-----------------------------------------
-
chandra_aca.maude_decom.
get_aca_images
(start, stop)[source]¶ Fetch ACA image telemetry
Parameters: - start – timestamp interpreted as a Chandra.Time.DateTime
- stop – timestamp interpreted as a Chandra.Time.DateTime
Returns: astropy.table.Table
-
chandra_aca.maude_decom.
get_aca_packets
(start, stop, level0=False, combine=False, adjust_time=False, calibrate=False)[source]¶ Fetch VCDU 1025-byte frames, extract ACA packets, unpack them and store them in a table.
Incomplete ACA packets (if there is a minor frame missing) can be combined or not into records with complete ACA telemetry. Compare these to calls to the function:
>>> from chandra_aca import maude_decom >>> img = maude_decom.get_aca_packets(684089000, 684089016, combine=True) >>> img = img[img['IMGNUM'] == 0] >>> img['TIME', 'MJF', 'MNF', 'COMMCNT', 'GLBSTAT', 'IMGTYPE', 'IMGROW0', 'IMGCOL0', >>> 'TEMPCCD', 'TEMPHOUS'] <Table masked=True length=4> TIME MJF MNF COMMCNT GLBSTAT IMGTYPE IMGROW0 IMGCOL0 TEMPCCD TEMPHOUS float64 uint32 uint32 uint8 uint8 uint8 int16 int16 int16 int16 ------------- ------ ------ ------- ------- ------- ------- ------- ------- -------- 684089001.869 78006 32 0 0 4 469 -332 -20 83 684089005.969 78006 48 0 0 4 469 -332 -20 83 684089010.069 78006 64 0 0 4 469 -332 -20 83 684089014.169 78006 80 0 0 4 469 -332 -20 83
Using combined=False, results in records with incomplete images. In this case, data can be missing from some records. For example, with 8X8 images, IMGROW0 and IMGCOL0 are present in the first ACA packet (image type 4) while the temperature is present in the second (image type 5):
>>> from chandra_aca import maude_decom >>> img = maude_decom.get_aca_packets(684089000, 684089016, combine=False) >>> img = img[img['IMGNUM'] == 0] >>> img['TIME', 'MJF', 'MNF', 'COMMCNT', 'GLBSTAT', 'IMGTYPE', 'IMGROW0', 'IMGCOL0', >>> 'TEMPCCD', 'TEMPHOUS'] <Table masked=True length=15> TIME MJF MNF COMMCNT GLBSTAT IMGTYPE IMGROW0 IMGCOL0 TEMPCCD TEMPHOUS float64 uint32 uint32 uint8 uint8 uint8 int16 int16 int16 int16 ------------- ------ ------ ------- ------- ------- ------- ------- ------- -------- 684089000.844 78006 28 0 0 7 -- -- -- -- 684089001.869 78006 32 0 0 4 469 -332 -- -- 684089002.894 78006 36 0 0 5 -- -- -20 83 684089003.919 78006 40 0 0 6 -- -- -- -- 684089004.944 78006 44 0 0 7 -- -- -- -- 684089005.969 78006 48 0 0 4 469 -332 -- -- 684089006.994 78006 52 0 0 5 -- -- -20 83 684089008.019 78006 56 0 0 6 -- -- -- -- 684089009.044 78006 60 0 0 7 -- -- -- -- 684089010.069 78006 64 0 0 4 469 -332 -- -- 684089011.094 78006 68 0 0 5 -- -- -20 83 684089012.119 78006 72 0 0 6 -- -- -- -- 684089013.144 78006 76 0 0 7 -- -- -- -- 684089014.169 78006 80 0 0 4 469 -332 -- -- 684089015.194 78006 84 0 0 5 -- -- -20 83
>>> img['IMG'].data[1] masked_BaseColumn(data = [[60.0 97.0 70.0 120.0 74.0 111.0 103.0 108.0] [67.0 90.0 144.0 96.0 88.0 306.0 82.0 67.0] [-- -- -- -- -- -- -- --] [-- -- -- -- -- -- -- --] [-- -- -- -- -- -- -- --] [-- -- -- -- -- -- -- --] [-- -- -- -- -- -- -- --] [-- -- -- -- -- -- -- --]], mask = [[False False False False False False False False] [False False False False False False False False] [ True True True True True True True True] [ True True True True True True True True] [ True True True True True True True True] [ True True True True True True True True] [ True True True True True True True True] [ True True True True True True True True]], fill_value = 1e+20)
>>> img['IMG'].data[2] masked_BaseColumn(data = [[-- -- -- -- -- -- -- --] [-- -- -- -- -- -- -- --] [76.0 81.0 160.0 486.0 449.0 215.0 88.0 156.0] [68.0 91.0 539.0 483.0 619.0 412.0 105.0 77.0] [-- -- -- -- -- -- -- --] [-- -- -- -- -- -- -- --] [-- -- -- -- -- -- -- --] [-- -- -- -- -- -- -- --]], mask = [[ True True True True True True True True] [ True True True True True True True True] [False False False False False False False False] [False False False False False False False False] [ True True True True True True True True] [ True True True True True True True True] [ True True True True True True True True] [ True True True True True True True True]], fill_value = 1e+20)
>>> img['IMG'].data[3] masked_BaseColumn(data = [[-- -- -- -- -- -- -- --] [-- -- -- -- -- -- -- --] [-- -- -- -- -- -- -- --] [-- -- -- -- -- -- -- --] [86.0 101.0 408.0 344.0 556.0 343.0 122.0 67.0] [196.0 195.0 114.0 321.0 386.0 115.0 69.0 189.0] [-- -- -- -- -- -- -- --] [-- -- -- -- -- -- -- --]], mask = [[ True True True True True True True True] [ True True True True True True True True] [ True True True True True True True True] [ True True True True True True True True] [False False False False False False False False] [False False False False False False False False] [ True True True True True True True True] [ True True True True True True True True]], fill_value = 1e+20)
>>> img['IMG'].data[4] Out[10]: masked_BaseColumn(data = [[-- -- -- -- -- -- -- --] [-- -- -- -- -- -- -- --] [-- -- -- -- -- -- -- --] [-- -- -- -- -- -- -- --] [-- -- -- -- -- -- -- --] [-- -- -- -- -- -- -- --] [67.0 61.0 67.0 176.0 99.0 72.0 79.0 88.0] [70.0 62.0 101.0 149.0 163.0 89.0 60.0 76.0]], mask = [[ True True True True True True True True] [ True True True True True True True True] [ True True True True True True True True] [ True True True True True True True True] [ True True True True True True True True] [ True True True True True True True True] [False False False False False False False False] [False False False False False False False False]], fill_value = 1e+20)
Parameters: - start – timestamp interpreted as a Chandra.Time.DateTime
- stop – timestamp interpreted as a Chandra.Time.DateTime
- level0 – bool. Implies combine=True, adjust_time=True, calibrate=True
- combine – bool. If True, multiple ACA packets are combined to form an image (depending on size), If False, ACA packets are not combined, resulting in multiple lines for 6x6 and 8x8 images.
- adjust_time – bool If True, half the integration time is subtracted
- calibrate – bool If True, pixel values will be ‘value * imgscale / 32 - 50’ and temperature values will be: 0.4 * value + 273.15
Returns: astropy.table.Table
-
chandra_aca.maude_decom.
get_raw_aca_packets
(start, stop)[source]¶ Fetch 1025-byte VCDU frames using maude and extract a list of 225-byte ACA packets.
If the first minor frame in a group of four ACA packets is within (start, stop), the three following minor frames are included if present.
returns a dictionary with keys [‘TIME’, ‘MNF’, ‘MJF’, ‘packets’, ‘flags’]. These correspond to the minor frame time, minor frame count, major frame count, the list of packets, and flags returned by Maude respectively.
Parameters: - start – timestamp interpreted as a Chandra.Time.DateTime
- stop – timestamp interpreted as a Chandra.Time.DateTime
Returns: