AstroIntensityMap

class sofia_redux.scan.source_models.astro_intensity_map.AstroIntensityMap(info, reduction=None)[source]

Bases: AstroData2D

Initialize an astronomical intensity map.

The astronomical intensity map represents the source as an Observation2D map containing data, noise, and exposure values. It also contains a base image containing the results of the previous reduction iteration in order to calculate gain and coupling increments.

Parameters:
infosofia_redux.scan.info.info.Info

The Info object which should belong to this source model.

reductionsofia_redux.scan.reduction.reduction.Reduction, optional

The reduction for which this source model should be applied.

Attributes Summary

flagspace

Return the flagspace associated with the underlying map observation.

referenced_attributes

Return attributes that should be referenced during a copy.

shape

Return the shape of the map.

Methods Summary

add_base()

Add the base to the observation map.

add_frames_from_integration(integration, ...)

Add frames from an integration to the source model.

add_model_data(source_model[, weight])

Add an increment source model data onto the current model.

add_points(frames, pixels, frame_gains, ...)

Add points to the source model.

base_footprint(pixels)

Returns the base footprint.

calculate_coupling(integration, pixels, ...)

Calculate the channel source coupling values.

clear_all_memory()

Clear all memory references prior to deletion.

copy([with_contents])

Return a copy of the source model.

covariant_points()

Return the number of points in the smoothing beam of the map.

create_from(scans[, assign_scans])

Initialize model from scans.

create_map()

Create the source model map.

filter_beam_correct()

Apply beam filter correction.

filter_source(filter_fwhm[, ...])

Apply source filtering.

get_clean_local_copy([full])

Get an unprocessed copy of the source model.

get_data()

Return the map data.

get_jansky_unit()

Return the Jansky's per unit area.

get_map_2d()

Return the 2D map.

get_peak_coords()

Return the coordinates of the peak value.

get_peak_index()

Return the peak index.

get_peak_source([degree, reduce_degrees])

Return the peak source model.

get_pixel_footprint()

Returns the Number of bytes per pixel.

get_source_name()

Return the source name.

get_table_entry(name)

Return a parameter value for a given name.

get_unit()

Return the map data unit.

is_masked()

Return the map mask.

mask_integration_samples(integration[, flag])

Set the sample flag in an integration for masked map elements.

mask_samples([flag])

Propagate masked source samples to integration sample flags.

mem_correct(lg_multiplier)

Apply maximum entropy correction to the map.

merge_accumulate(other)

Merge another source with this one.

merge_mask(other_map)

Merge the mask from another map onto this one.

post_process_scan(scan)

Perform post-processing steps on a scan.

process_final()

Perform the final processing steps.

reset_filtering()

Reset the map filtering.

reset_processing()

Reset the source processing.

set_base()

Set the base to the map (copy of).

set_data_shape(shape)

Set the shape of the map.

set_filtering(fwhm)

Set the map filtering FWHM.

set_info(info)

Set the Info object for the source model.

smooth_to(fwhm)

Smooth the map to a given FWHM.

stand_alone()

Create a stand alone base image.

sync_source_gains(frames, pixels, ...)

Remove the map source from frame data.

update_mask([blanking_level, min_neighbors])

Update the map mask based on significance levels and valid neighbors.

Attributes Documentation

flagspace

Return the flagspace associated with the underlying map observation.

Returns:
sofia_redux.scan.flags.flags.Flag
referenced_attributes

Return attributes that should be referenced during a copy.

Returns:
set (str)
shape

Return the shape of the map.

Note that this is in numpy (y, x) order.

Returns:
tuple of int

Methods Documentation

add_base()[source]

Add the base to the observation map.

Returns:
None
add_frames_from_integration(integration, pixels, source_gains, signal_mode=None)[source]

Add frames from an integration to the source model.

Parameters:
integrationIntegration

The integration to add.

pixelsChannelGroup

The channels (pixels) to add to the source model.

source_gainsnumpy.ndarray (float)

The source gains for the all channels (pixels). Should be of shape (all_channels,).

signal_modeFrameFlagTypes, optional

The signal mode flag, indicating which signal should be used to extract the frame source gains. Typically, TOTAL_POWER.

Returns:
mapping_framesint

The number of frames that contributed towards mapping.

add_model_data(source_model, weight=1.0)[source]

Add an increment source model data onto the current model.

Parameters:
source_modelAstroIntensityMap

The source model increment.

weightfloat, optional

The weight of the source model increment.

Returns:
None
add_points(frames, pixels, frame_gains, source_gains)[source]

Add points to the source model.

Accumulates the combined frame and channel data to the source map for each frame/channel sample. If a given sample maps to a single map pixel, that pixel is incremented by:

i = frame_data * weights * gains
w = weights * gains^2
weights = frame_weight / channel_variance
gains = frame_gain * channel_gain

where i is the weighted data increment, and w is the weight increment. The exposure values are also added to by simply incrementing the time at any pixel by the sampling interval multiplied by the number of samples falling within that pixel.

Parameters:
framesFrames

The frames to add to the source model.

pixelsChannelGroup

The channels (pixels) to add to the source model.

frame_gainsnumpy.ndarray (float)

The gain values for all frames of shape (n_frames,).

source_gainsnumpy.ndarray (float)

The channel source gains for all channels of shape (all_channels,).

Returns:
mapping_framesint

The number of valid mapping frames added for the model.

base_footprint(pixels)[source]

Returns the base footprint.

Parameters:
pixelsint
Returns:
int
calculate_coupling(integration, pixels, source_gains, sync_gains)[source]

Calculate the channel source coupling values.

Parameters:
integrationIntegration
pixelsChannels or ChannelGroup

The pixels for which to calculate coupling.

source_gainsnumpy.ndarray (float)

The source gains for all frames in the integration of shape (n_frames,).

sync_gainsnumpy.ndarray (float)

The sync gains for all channels in the integration of shape (all_channels,).

Returns:
None
clear_all_memory()[source]

Clear all memory references prior to deletion.

Returns:
None
copy(with_contents=True)[source]

Return a copy of the source model.

Parameters:
with_contentsbool, optional

If True, return a true copy of the map. Otherwise, just return a map with basic metadata.

Returns:
AstroIntensityMap
covariant_points()[source]

Return the number of points in the smoothing beam of the map.

Returns:
float
create_from(scans, assign_scans=True)[source]

Initialize model from scans.

Sets the model scans to those provided, and the source model for each scan as this. All integration gains are normalized to the first scan. If the first scan is non-sidereal, the system will be forced to an equatorial frame.

Parameters:
scanslist (Scan)

A list of scans from which to create the model.

assign_scansbool, optional

If True, assign the scans to this source model. Otherwise, there will be no hard link between the scans and source model.

Returns:
None
create_map()[source]

Create the source model map.

Returns:
None
filter_beam_correct()[source]

Apply beam filter correction.

Returns:
None
filter_source(filter_fwhm, filter_blanking=None, use_fft=False)[source]

Apply source filtering.

Parameters:
filter_fwhmastropy.units.Quantity or Coordinate

The filtering FWHM for the source.

filter_blankingfloat, optional

Only apply filtering within the optional range -filter_blanking -> +filter_blanking.

use_fftbool, optional

If True, use FFTs to perform the filtering. Otherwise, do full convolutions.

Returns:
None
get_clean_local_copy(full=False)[source]

Get an unprocessed copy of the source model.

Parameters:
fullbool, optional

If True, copy additional parameters for stand-alone reductions that would otherwise be referenced.

Returns:
SourceModel
get_data()[source]

Return the map data.

Returns:
Observation2D
get_jansky_unit()[source]

Return the Jansky’s per unit area.

Returns:
janskysunits.Quantity
get_map_2d()[source]

Return the 2D map.

Returns:
Map2D
get_peak_coords()[source]

Return the coordinates of the peak value.

Returns:
peak_coordinatesCoordinate2D

The (x, y) peak coordinate.

get_peak_index()[source]

Return the peak index.

Returns:
indexIndex2D

The peak (x, y) coordinate.

get_peak_source(degree=3, reduce_degrees=False)[source]

Return the peak source model.

Parameters:
degreeint, optional

The spline degree used to fit the peak map value.

reduce_degreesbool, optional

If True, allow the spline fit to reduce the number of degrees in cases where there are not enough points available to perform the spline fit of degree. If False, a ValueError will be raised if such a fit fails.

Returns:
GaussianSource
get_pixel_footprint()[source]

Returns the Number of bytes per pixel.

This is probably no longer relevant and just copies CRUSH.

Returns:
int
get_source_name()[source]

Return the source name.

Returns:
str
get_table_entry(name)[source]

Return a parameter value for a given name.

Parameters:
namestr, optional
Returns:
value
get_unit()[source]

Return the map data unit.

Returns:
astropy.units.Quantity
is_masked()[source]

Return the map mask.

Returns:
numpy.ndarray (bool)
mask_integration_samples(integration, flag='SAMPLE_SKIP')[source]

Set the sample flag in an integration for masked map elements.

Parameters:
integrationIntegration
flagstr or int or Enum, optional

The name, integer identifier, or actual flag by which to flag samples.

Returns:
None
mask_samples(flag='SAMPLE_SKIP')[source]

Propagate masked source samples to integration sample flags.

Flags the integration samples using the given flag if they correspond to a masked source sample.

Parameters:
flagstr or int or Enum, optional

The name, integer identifier, or actual flag by which to mask samples. By default, any integration samples that correspond to masked source samples will be flagged as ‘SAMPLE_SKIP’.

Returns:
None
mem_correct(lg_multiplier)[source]

Apply maximum entropy correction to the map.

Parameters:
lg_multiplierfloat
Returns:
None
merge_accumulate(other)[source]

Merge another source with this one.

Parameters:
otherAstroIntensityMap
Returns:
None
merge_mask(other_map)[source]

Merge the mask from another map onto this one.

Parameters:
other_mapFlaggedArray
Returns:
None
post_process_scan(scan)[source]

Perform post-processing steps on a scan.

At this stage, the map should have already been added to the main reduction map and we are left with a throw-away map that can be used to update settings in other objects. The intensity map may be used to update the pointing for a given scan.

Parameters:
scanScan
Returns:
None
process_final()[source]

Perform the final processing steps.

The additional steps performed for the AstroIntensityMap are map leveling (if not extended or deep) and map re-weighting. The map may also be resampled if re-griding is enabled.

Returns:
None
reset_filtering()[source]

Reset the map filtering.

Returns:
None
reset_processing()[source]

Reset the source processing.

Returns:
None
set_base()[source]

Set the base to the map (copy of).

Returns:
None
set_data_shape(shape)[source]

Set the shape of the map.

Parameters:
shapetuple (int)
Returns:
None
set_filtering(fwhm)[source]

Set the map filtering FWHM.

Parameters:
fwhmastropy.units.Quantity
Returns:
None
set_info(info)[source]

Set the Info object for the source model.

This sets the provided info as the primary Info object containing the configuration and reduction information for the source model. The source model will also take ownership of the info and set various parameters from the contents.

Parameters:
infosofia_redux.info.info.Info
Returns:
None
smooth_to(fwhm)[source]

Smooth the map to a given FWHM.

Parameters:
fwhmastropy.units.Quantity or Gaussian2D
Returns:
None
stand_alone()[source]

Create a stand alone base image.

Returns:
None
sync_source_gains(frames, pixels, frame_gains, source_gains, sync_gains)[source]

Remove the map source from frame data.

In addition to source removal, samples are also flagged if the map is masked at that location.

For a given sample at frame i and channel j, frame data d_{i,j} will be decremented by dg where:

dg = fg * ( (gain(source) * map[index]) - (gain(sync) * base[index]) )

Here, fg is the frame gain and index is the index on the map of sample (i,j).

Any masked map value will result in matching samples being flagged.

Parameters:
framesFrames

The frames for which to remove the source gains.

pixelsChannels or ChannelGroup

The channels for which to remove the source gains.

frame_gainsnumpy.ndarray (float)

An array of frame gains of shape (n_frames,).

source_gainsnumpy.ndarray (float)

An array of channel source gains for all channels of shape (all_channels,).

sync_gainsnumpy.ndarray (float)

an array of channel sync gains for all channels of shape (all_channels,). The sync gains should contain the prior source gain.

Returns:
None
update_mask(blanking_level=None, min_neighbors=2)[source]

Update the map mask based on significance levels and valid neighbors.

If a blanking level is supplied, significance values above or equal to the blanking level will be masked.

Parameters:
blanking_levelfloat, optional

The significance level used to mark the map. If not supplied, significance is irrelevant. See above for more details.

min_neighborsint, optional

The minimum number of neighbors including the pixel itself. Therefore, the default of 2 excludes single pixels as this would require a single valid pixel and one valid neighbor.

Returns:
None