corgidrp.l3_to_l4#
Functions#
|
Interpolate over bad pixels in image and error arrays using a median filter. |
|
Applies the distortion correction to the dataset. The function assumes bad |
|
Determines the star position within a coronagraphic dataset by analyzing frames that |
|
Perform PSF subtraction on the dataset. Optionally using a reference star dataset. |
|
Derotate the Image, ERR, and DQ data by the angle offset to make the FoV up to North. |
|
A procedure for estimating the centroid of the zero-point image |
|
add_wavelength_map adds the wavelength map + error and the position lookup table as extensions to the frames |
|
Find the position of the star using the information from the satellite |
|
extract an optionally error weighted 1D - spectrum and wavelength information of a point source from a box around |
|
Aligns a dataset of 2D images by recentering them using the STARLOCX and STARLOCY header keywords |
|
Aligns the frames by centering them on STARLOC |
|
Takes in polarimetric L3 images and their unocculted polarimetric observations, |
|
Takes in a L3 polarimetric dataset, performs PSF subtraction on total intensity, calculates |
|
RDI PSF subtraction for spectroscopy mode. |
|
Calculates spectroscopy PSF subtraction algorithmic throughput. |
|
combination of psf subtracted frames for spectroscopy mode. |
|
Updates the data level to L4. Only works on L3 data. |
|
Updates the data level to L4 for polarimetry data. Only works on L3 data. |
Module Contents#
- corgidrp.l3_to_l4.replace_bad_pixels(input_dataset, kernelsize=3, dq_thresh=1)[source]#
Interpolate over bad pixels in image and error arrays using a median filter. TODO: Add additional options for bad pixel replacement (e.g. 2d interpolation that can handle nans, constant value, etc.)
- Parameters:
input_dataset (corgidrp.data.Dataset) – input L3 dataset with bad pixels
kernelsize (int, optional) – Size of median filter window in pixels. Defaults to 3.
dq_thresh (int, optional) – Minimum DQ value for a pixel to be replaced. Defaults to 1.
- Returns:
A copy of the dataset with the bad pixels and error values interpolated over. The bad pixel map is unchanged.
- Return type:
- corgidrp.l3_to_l4.distortion_correction(input_dataset, astrom_calibration)[source]#
Applies the distortion correction to the dataset. The function assumes bad pixels have been corrected beforehand. Furthermore it also applies the distortion correction to the error maps.
- Parameters:
input_dataset (corgidrp.data.Dataset) – a dataset of Images (L3-level)
astrom_calibration (corgidrp.data.AstrometricCalibration) – an AstrometricCalibration calibration file to model the distortion
- Returns:
a version of the input dataset with the distortion correction applied
- Return type:
- corgidrp.l3_to_l4.find_star(input_dataset, star_coordinate_guess=None, thetaOffsetGuess=0, satellite_spot_parameters=None, drop_satspots_frames=True, pri_split_keywords=None, ext_split_keywords=None)[source]#
Determines the star position within a coronagraphic dataset by analyzing frames that contain satellite spots (indicated by
SATSPOTS=Truein the image header). The function computes the median of all science frames (SATSPOTS=False) and the median of all satellite spot frames (SATSPOTS=True), then estimates the star location based on these median images and the initial guess provided.The star’s (x, y) location is stored in each frame’s extension header under
STARLOCXandSTARLOCY.In case of polarimetric data, the star location is estimated on the first slice and the second slice is aligned on it. POL 0 and POL 45 are processed independantly
You can replace many of the default settings for by adjusting the satellite_spot_parameters dictionary. You only need to replace the parameters of interest and the rest will stay as defaults.
- satellite_spot_parameters of the form:
- offsetdict
Parameters for estimating the offset of the star center:
- spotSepPixfloat
Expected (model-based) separation of the satellite spots from the star. Units: pixels.
- roiRadiusPixfloat
Radius of the region of interest around each satellite spot. Units: pixels.
- probeRotVecDegarray_like
Angles (degrees CCW from x-axis) specifying the position of satellite spot pairs.
- nSubpixelsint
Number of subpixels across for calculating region-of-interest mask edges.
- nStepsint
Number of points in grid search along each direction.
- stepSizefloat
Step size for the grid search. Units: pixels.
- nIterint
Number of iterations refining the radial separation.
- separationdict
Parameters for estimating the separation of satellite spots from the star:
- spotSepPixfloat
Expected separation between star and satellite spots. Units: pixels.
- roiRadiusPixfloat
Radius of the region of interest around each satellite spot. Units: pixels.
- probeRotVecDegarray_like
Angles (degrees CCW from x-axis) specifying the position of satellite spot pairs.
- nSubpixelsint
Number of subpixels across for calculating region-of-interest mask edges.
- nStepsint
Number of points in grid search along each direction.
- stepSizefloat
Step size for the grid search. Units: pixels.
- nIterint
Number of iterations refining the radial separation.
- Parameters:
input_dataset (corgidrp.data.Dataset) – A dataset of L3-level frames. Frames should be labeled in their primary headers with
SATSPOTS=False(science frames) orSATSPOTS=True(satellite spot frames).star_coordinate_guess (tuple of float or None, optional) – Initial guess for the star’s (x, y) location as absolute coordinates. If
None, defaults to the center of the median satellite spot image. Defaults to None.thetaOffsetGuess (float, optional) – Initial guess for any angular rotation of the star center (in degrees, for example). Defaults to 0.
satellite_spot_parameters (dict, optional) – Dictionary containing tuning parameters for spot separation and offset estimation. The dictionary can contain the following keys and structure. Only provided parameters will be changed, otherwise defaults for the mode will be used: If None, default parameters corresponding to the specified observing_mode will be used.
drop_satspots_frames (bool, optional) – If True, frames with satellite spots (
SATSPOTS=True) will be removed from the returned dataset. Defaults to True.pri_split_keywords (list of str, optional) – List of primary header keywords to use for splitting the dataset into subsets. If None, defaults to [‘VISITID’]. Defaults to None.
ext_split_keywords (list of str, optional) – List of extension header keywords to use for splitting the dataset into subsets. If None, defaults to [‘DPAMNAME’]. Defaults to None.
- Returns:
The original dataset, augmented with the star’s (x, y) location stored in the extension header (
ext_hdr) of each frame under the keysSTARLOCXandSTARLOCY.- Return type:
- Raises:
AssertionError – If any frames have an invalid
SATSPOTSkeyword (not 0 or 1), or if the frames do not all share the same observing mode (as determined by theFSMPRFLkeyword).
Notes
This function merges the science frames (for reference) and the satellite spot frames (for analysis) by taking a median image of each set.
The star center is computed using the median images and the
star_center.star_center_from_satellite_spotsroutine.Future enhancements may include separate handling of positive vs. negative satellite spot frames once the relevant metadata keywords are defined.
This routine can fail, if the guess position is off by more than a few pixels. More than 2 pixels on any axis leads almost systematically to failure A significantly wrong guess of the angle offset can also lead to failure.
- corgidrp.l3_to_l4.do_psf_subtraction(input_dataset, ct_calibration=None, reference_star_dataset=None, outdir=None, fileprefix='', measure_klip_thrupt=True, measure_1d_core_thrupt=True, cand_locs=None, kt_seps=None, kt_pas=None, kt_snr=20.0, num_processes=None, dq_thresh=1, **klip_kwargs)[source]#
Perform PSF subtraction on the dataset. Optionally using a reference star dataset. .. todo:
Propagate error correctly Use corgidrp.combine.combine_images() to do time collapse? What info is missing from output dataset headers? Add comments to new ext header cards Make sure psfsub test output data gets saved in a reasonable place Update output filename to: CGI_<Last science target VisitID>_<Last science target TimeUTC>_L<>.fits
- Parameters:
input_dataset (corgidrp.data.Dataset) – a dataset of Images (L3-level)
ct_calibration (corgidrp.data.CoreThroughputCalibration, optional) – core throughput calibration object. Required if measuring KLIP throughput or 1D core throughput. Defaults to None.
reference_star_dataset (corgidrp.data.Dataset, optional) – a dataset of Images of the reference star. If not provided, references will be searched for in the input dataset.
outdir (str or path, optional) – path to output directory. Defaults to “KLIP_SUB”.
fileprefix (str, optional) – prefix of saved output files. Defaults to “”.
measure_klip_thrupt (bool, optional) – Whether to measure KLIP throughput via injection-recovery. Separations and throughput levels for each separation and KL mode are saved in Dataset[0].hdu_list[‘KL_THRU’]. Defaults to True.
measure_1d_core_thrupt (bool, optional) – Whether to measure the core throughput as a function of separation. Separations and throughput levels for each separation are saved in Dataset[0].hdu_list[‘CT_THRU’]. Defaults to True.
cand_locs (list of tuples, optional) – Locations of known off-axis sources, so we don’t inject a fake PSF too close to them. This is a list of tuples (sep_pix,pa_degrees) for each source. Defaults to [].
kt_seps (np.array, optional) – Separations (in pixels from the star center) at which to inject fake PSFs for KLIP throughput calibration. If not provided, a linear spacing of separations between the IWA & OWA will be chosen.
kt_pas (np.array, optional) – Position angles (in degrees counterclockwise from north/up) at which to inject fake PSFs at each separation for KLIP throughput calibration. Defaults to [0.,90.,180.,270.].
kt_snr (float, optional) – SNR of fake signals to inject during KLIP throughput calibration. Defaults to 20.
num_processes (int) – number of processes for parallelizing the PSF subtraction
dq_thresh (int) – DQ threshold for considering a pixel bad. Defaults to 1.
klip_kwargs – Additional keyword arguments to be passed to pyKLIP fm.klip_dataset, as defined here <https://pyklip.readthedocs.io/en/latest/pyklip.html#pyklip.fm.klip_dataset>. ‘mode’, e.g. ADI/RDI/ADI+RDI, is chosen autonomously if not specified. ‘annuli’ defaults to 1. ‘annuli_spacing’ defaults to ‘constant’. ‘subsections’ defaults to 1. ‘movement’ defaults to 1. ‘numbasis’ defaults to [1,4,8,16], ‘time_collapse’ defaults to ‘mean’.
- Returns:
a version of the input dataset with the PSF subtraction applied (L4-level)
- Return type:
- corgidrp.l3_to_l4.northup(input_dataset, use_wcs=True, rot_center='im_center', new_center=None)[source]#
Derotate the Image, ERR, and DQ data by the angle offset to make the FoV up to North. The northup function looks for ‘STARLOCX’ and ‘STARLOCY’ for the star location. If not, it uses the center of the FoV as the star location. With use_wcs=True it uses WCS infomation to calculate the north position angle, or use just ‘PA_APER’ header keyword if use_wcs is False (not recommended). TODO: Update pixel locations that are saved in the header! TODO: Add tests for behavior of new_center
- Parameters:
input_dataset (corgidrp.data.Dataset) – a dataset of Images (L3-level) - now handles pol datasets shapes
use_wcs – if you want to use WCS to correct the north position angle, set True (default). rot_center: ‘im_center’, ‘starloc’, or manual coordinate (x,y). ‘im_center’ uses the center of the image. ‘starloc’ refers to ‘STARLOCX’ and ‘STARLOCY’ in the header.
new_center – location (xy) to move the center to after rotation.
- Returns:
North is up, East is left
- Return type:
- corgidrp.l3_to_l4.determine_wave_zeropoint(input_dataset, spec_filter_offset, template_dataset=None, xcent_guess=None, ycent_guess=None, bb_nb_dx=None, bb_nb_dy=None, return_all=False)[source]#
A procedure for estimating the centroid of the zero-point image (satellite spot or PSF) taken through the narrowband filter (2C or 3D) and slit.
- Parameters:
input_dataset (corgidrp.data.Dataset) – Dataset containing 2D PSF or satellite spot images taken through the narrowband filter and slit.
spec_filter_offset (corgidrp.data.SpecFilterOffset) – instance of SpecFilterOffset calibration class
template_dataset (corgidrp.data.Dataset) – dataset of the template PSF, if None, a simulated PSF from the data/spectroscopy/template path is taken
xcent_guess (float) – initial x guess for the centroid fit for all frames
ycent_guess (float) – initial y guess for the centroid fit for all frames
bb_nb_dx (float) – horizontal image offset between the narrowband and broadband filters, in EXCAM pixels. This will override the offset in the existing lookup table.
bb_nb_dy (float) – vertical image offset between the narrowband and broadband filters, in EXCAM pixels. This will override the offset in the existing lookup table.
return_all (boolean) – if false (default) returns only the broad band science frames, if true it returns all (including narrow band) frames
- Returns:
- the returned science dataset without the satellite spots images and the wavelength zeropoint
information as header keywords, which is WAVLEN0, WV0_X, WV0_XERR, WV0_Y, WV0_YERR, WV0_DIMX, WV0_DIMY
- Return type:
- corgidrp.l3_to_l4.add_wavelength_map(input_dataset, disp_model, pixel_pitch_um=13.0, ntrials=1000)[source]#
add_wavelength_map adds the wavelength map + error and the position lookup table as extensions to the frames
- Parameters:
input_dataset (corgidrp.data.Dataset) – a dataset of spectroscopy Images (L3-level)
disp_model (corgidrp.data.DispersionModel) – dispersion model of the corresponding band
pixel_pitch_um (float) – EXCAM pixel pitch in microns, default: 13.0
ntrials (int) – number of trials when applying a Monte Carlo error propagation to estimate the uncertainties of the values in the wavelength calibration map
- Returns:
dataset with appended wavelength map and error
- Return type:
- corgidrp.l3_to_l4.find_spec_star(input_dataset, r_lamD=3, phi_deg=0)[source]#
Find the position of the star using the information from the satellite spot. The position of the satellite spot on EXCAM is given by the zero-point solution. Using the information of the commanded position of the satellite spot with respect the occulted star, one can infer the location of the occulted star. The relative of the satellite spot with respect the occulted star is given in polar coordinates. The radial distance of the satellite spot is measured in units lambda/D, with lambda the band reference wavelength, either 730 nm (band 3) or 660 nm (band 2), and D=2.4 m. The polar angle is measured in degrees, with 0 degrees meaning +X and 90 degrees meaning +Y. The polar coordinates of the satellite spot are translated into (X,Y) EXCAM pixel coordinates, which can then be subtracted from the zero-point solution to infer the location of the occulted star.
- Parameters:
input_dataset (corgidrp.data.Dataset) – a dataset of spectroscopy Images (L3-level)
r_lamD (float) – Radial distance of the satellite spot on EXCAM with respect
lambda/D. (the occulted star in units of)
phi_deg (float) – Polar angle of the satellite spot on EXCAM with respect
degrees (with 0 degrees meaning +X and 90)
degrees
+Y. (meaning)
- Returns:
- Dataset with updated keywords recording the satellite position
in EXCAM pixels.
- Return type:
- corgidrp.l3_to_l4.extract_spec(input_dataset, halfwidth=2, halfheight=9, apply_weights=False)[source]#
extract an optionally error weighted 1D - spectrum and wavelength information of a point source from a box around the wavelength zero point with units photoelectron/s/bin.
- Parameters:
input_dataset (corgidrp.data.Dataset)
halfwidth (int) – The width of the fitting region is 2 * halfwidth + 1 pixels across dispersion
halfheight (int) – The height of the fitting region is 2 * halfheight + 1 pixels along dispersion.
apply_weights (boolean) – if true a weighted sum is calculated using 1/error^2 as weights.
- Returns:
dataset containing the spectral 1D data, error and corresponding wavelengths
- Return type:
- corgidrp.l3_to_l4.align_2d_frames(input_dataset, center='first_frame')[source]#
Aligns a dataset of 2D images by recentering them using the STARLOCX and STARLOCY header keywords
- Parameters:
input_dataset (corgidrp.data.Dataset) – the L3-level dataset of 2D images with STARLOCX and STARLOCY
center (str or tuple) – Can be one of three options. 1. ‘first_frame’ (default) - aligns all frames to the STARLOCX/Y of the first frame 2. ‘im_center’ - aligns all frames to the center of the image 3. (x,y) tuple - aligns all frames to the provided (x,y) pixel location
- Returns:
L3 dataset where all the images are registered to the same pixel
- Return type:
- corgidrp.l3_to_l4.align_polarimetry_frames(input_dataset)[source]#
Aligns the frames by centering them on STARLOC
- Parameters:
input_dataset (corgidrp.data.Dataset) – the L3-level dataset of polarimetry images with STARLOCX and STARLOCY
- Returns:
L3 dataset where all the images are registered to the same pixel
- Return type:
- corgidrp.l3_to_l4.subtract_stellar_polarization(input_dataset, system_mueller_matrix_cal, nd_mueller_matrix_cal)[source]#
Takes in polarimetric L3 images and their unocculted polarimetric observations, computes and subtracts off the stellar polarization component from each image TODO: make issue about error propagation, need to check that it is done correctly
and make changes if necessary to ensure the errors are accurate
- Parameters:
input_dataset (corgidrp.data.Dataset) – a dataset of L3 images, must include unocculted observations taken with both wollastons at the same telescope rotation angle. All frames for the same target star must have the same x and y dimensions
system_mueller_matrix_cal (corgidrp.data.MuellerMatrix) – mueller matrix calibration of the system without a ND filter
nd_mueller_matrix_cal (corgidrp.data.NDMuellerMatrix) – mueller matrix calibration of the system with the ND filter used for unocculted observations
- Returns:
The input data with stellar polarization removed, excluding the unocculted observations
- Return type:
- corgidrp.l3_to_l4.combine_polarization_states(input_dataset, system_mueller_matrix_cal, ct_calibration, svd_threshold=1e-05, use_wcs=True, rot_center='im_center', reference_star_dataset=None, measure_klip_thrupt=True, measure_1d_core_thrupt=True, cand_locs=None, kt_seps=None, kt_pas=None, kt_snr=20.0, num_processes=None, **klip_kwargs)[source]#
Takes in a L3 polarimetric dataset, performs PSF subtraction on total intensity, calculates and returns the on-sky Stokes datacube TODO: determine how to propagate DQ extension through matrix inversion
- Parameters:
input_dataset (corgidrp.data.Dataset) – a dataset of polarimetric Images (L3-level), should be of the same size and same target
system_mueller_matrix_cal (corgidrp.data.MuellerMatrix) – mueller matrix calibration of the instrument
svd_threshold (float, optional) – The threshold for singular values in the SVD inversion. Defaults to 1e-5 (semi-arbitrary).
use_wcs (bool, optional) – Uses WCS coordinates to rotate northup, defaults to true. If false, uses PA_APER header instead.
rot_center (string, optional) – Define the center to rotate the images with respect to. Options are ‘im_center’, ‘starloc’, or manual coordinate (x,y). ‘im_center’ uses the center of the image. ‘starloc’ refers to ‘STARLOCX’ and ‘STARLOCY’ in the header.
ct_calibration (corgidrp.data.CoreThroughputCalibration, optional) – For PSF Subtraction. core throughput calibration object. Required if measuring KLIP throughput or 1D core throughput. Defaults to None.
reference_star_dataset (corgidrp.data.Dataset, optional) – For PSF Subtraction. a dataset of Images of the reference star. If not provided, references will be searched for in the input dataset.
measure_klip_thrupt (bool, optional) – For PSF Subtraction. Whether to measure KLIP throughput via injection-recovery. Separations and throughput levels for each separation and KL mode are saved in Dataset[0].hdu_list[‘KL_THRU’]. Defaults to False.
measure_1d_core_thrupt (bool, optional) – For PSF Subtraction. Whether to measure the core throughput as a function of separation. Separations and throughput levels for each separation are saved in Dataset[0].hdu_list[‘CT_THRU’]. Defaults to True.
cand_locs (list of tuples, optional) – For PSF Subtraction. Locations of known off-axis sources, so we don’t inject a fake PSF too close to them. This is a list of tuples (sep_pix,pa_degrees) for each source. Defaults to [].
kt_seps (np.array, optional) – For PSF Subtraction. Separations (in pixels from the star center) at which to inject fake PSFs for KLIP throughput calibration. If not provided, a linear spacing of separations between the IWA & OWA will be chosen.
kt_pas (np.array, optional) – For PSF Subtraction. Position angles (in degrees counterclockwise from north/up) at which to inject fake PSFs at each separation for KLIP throughput calibration. Defaults to [0.,90.,180.,270.].
kt_snr (float, optional) – For PSF Subtraction. SNR of fake signals to inject during KLIP throughput calibration. Defaults to 20.
num_processes (int) – For PSF Subtraction. number of processes for parallelizing the PSF subtraction
klip_kwargs – For PSF Subtraction. Additional keyword arguments to be passed to pyKLIP fm.klip_dataset, as defined here <https://pyklip.readthedocs.io/en/latest/pyklip.html#pyklip.fm.klip_dataset>. ‘mode’, e.g. ADI/RDI/ADI+RDI, is chosen autonomously if not specified. ‘annuli’ defaults to 1. ‘annuli_spacing’ defaults to ‘constant’. ‘subsections’ defaults to 1. ‘movement’ defaults to 1. ‘numbasis’ defaults to [1] if no reference star dataset is provided, and [len(reference_star_dataset)] otherwise, numbasis must be of length 1 only for this step function.
- Returns:
- Dataset consisting of a 4xnxm on-sky Stokes datacube, where the first dimension corresponds to IQUV
with I being the PSF-subtracted total intensity.
- Return type:
- corgidrp.l3_to_l4.spec_psf_subtraction(input_dataset, mask_pl_pixels=3, bg_threshold=5, poly_order=5, sigma_injplanet=1.5)[source]#
RDI PSF subtraction for spectroscopy mode. Assumes the reference images are marked with PSFREF=True in the primary header and that they all have the same alignment.
- Parameters:
input_dataset (corgidrp.data.Dataset) – L3 dataset containing the science and reference images.
mask_pl_pixels (int, optional) – No. of pixels to mask on either side of planet position. Default value is 3 pixels.
bg_threshold (int, optional) – Threshold to identify detector rows with reference psf using no. of std deviations from median. Default value is 5.
poly_order (int, optional) – Order of polynomial used to scale down reference star psf to target star psf. Default value is 5.
sigma_injplanet (float, optional) – Standard deviation of injected gaussian planet for throughput calcs from injection-recovery tests. Default value is 1.5 pixels.
- Returns:
dataset containing the PSF-subtracted science images
- Return type:
- corgidrp.l3_to_l4.spec_throughput(orig_frame, shifted_ref, injection_site, sigma_injplanet, best_peak_col, start_row, end_row, poly_order=5, amplitude=0.02, mask_pl_pixels=3)[source]#
Calculates spectroscopy PSF subtraction algorithmic throughput.
- Parameters:
orig_frame (corgidrp.data.Image) – Science data frame from spec_psf_subtraction.
shifted_ref (corgidrp.data.Image) – Shifted reference frame from spec_psf_subtraction.
injection_site (int) – Column with peak of injected planet psf.
sigma_injplanet (float) – Standard deviation of injected gaussian planet for throughput calcs from injection-recovery tests.
best_peak_col (int) – Column with peak of reference star speckles.
start_row (int) – Start row location of ref star speckles.
end_row (int) – End row location of ref star speckles.
poly_order (int, optional) – Order of polynomial used to scale down reference star psf to target star psf. Default value is 5.
amplitude (float, optional) – Amplitude of planet psf. Default value is 0.02.
mask_pl_pixels (int, optional) – No. of pixels to mask on either side of planet position. Also used to define pixels where planet is injected. Default value is 3 pixels.
- Returns:
Algorithmic throughput as a function of wavelength.
- corgidrp.l3_to_l4.combine_spec(input_dataset, collapse='mean', num_frames_scaling=True)[source]#
combination of psf subtracted frames for spectroscopy mode. Assumes they all have the same alignment.
- Parameters:
input_dataset (corgidrp.data.Dataset) – L3 spectroscopy dataset.
collapse (str) – “mean” or “median”. (default: mean)
num_frames_scaling (bool) – Multiply by number of frames in sequence in order to ~conserve photons (default: True)
- Returns:
dataset containing one frame of the PSF-subtracted science images
- Return type:
- corgidrp.l3_to_l4.update_to_l4(input_dataset, corethroughput_cal, flux_cal)[source]#
Updates the data level to L4. Only works on L3 data.
Currently only checks that data is at the L3 level
- Parameters:
input_dataset (corgidrp.data.Dataset) – a dataset of Images (L3-level)
corethroughput_cal (corgidrp.data.CoreThroughputCalibration) – a CoreThroughputCalibration calibration file. Can be None
flux_cal (corgidrp.data.FluxCalibration) – a FluxCalibration calibration file. Cannot be None
- Returns:
same dataset now at L4-level
- Return type:
- corgidrp.l3_to_l4.update_to_l4_pol(input_dataset, corethroughput_cal, flux_cal_pol0, flux_cal_pol45)[source]#
Updates the data level to L4 for polarimetry data. Only works on L3 data.
Calls update_to_l4 for shared logic, then replaces the single FLXCALFN header with separate FLXCLF0 (POL0) and FLXCLF45 (POL45) headers.
- Parameters:
input_dataset (corgidrp.data.Dataset) – a dataset of Images (L3-level)
corethroughput_cal (corgidrp.data.CoreThroughputCalibration) – a CoreThroughputCalibration calibration file. Can be None
flux_cal_pol0 (corgidrp.data.FluxcalFactor) – the POL0 FluxcalFactor calibration file
flux_cal_pol45 (corgidrp.data.FluxcalFactor) – the POL45 FluxcalFactor calibration file
- Returns:
same dataset now at L4-level
- Return type: