pyklip.fmlib package
Submodules
pyklip.fmlib.diskfm module
- class pyklip.fmlib.diskfm.DiskFM(inputs_shape, numbasis, dataset, model_disk, basis_filename='', kl_basis_file=None, load_from_basis=False, save_basis=False, aligned_center=None, psf_library=None, mode=None, annuli=None, subsections=None, numthreads=None)[source]
Bases:
NoFM
Defining a model disk to which we apply the Forward Modelling. There are 3 ways:
“Save Basis mode” (save_basis=true), we are preparing to save the FM basis
“Load Basis mode” (load_from_basis = true), most of the parameters are derived from the previous fm.klip_dataset which measured FM basis.
“Simple FM mode” (save_basis = load_from_basis = False). Just for a unique disk FM.
- Parameters:
inputs_shape – shape of the inputs numpy array. Typically (N, x, y)
numbasis – 1d numpy array consisting of the number of basis vectors to use
dataset – an instance of Instrument.Data. We need it to know the parameters to “prepare” first inital model.
model_disk – a model of the disk of size (wvs, x, y) or (x, y)
basis_filename – filename to save and load the KL basis. Filenames can haves 2 recognizable extensions: .h5 or .pkl. We strongly recommand .h5 as pickle have problem of compatibility between python 2 and 3 and sometimes between computer (e.g. KL modes not readable on another computer)
load_from_basis – if True, load the KL basis at basis_filename. It only need to be done once, after which you can measure FM with only update_model()
save_basis – if True, save the KL basis at basis_filename. If load_from_basis is True, save_basis is automatically set to False, it is useless to load and save the matrix at the same time.
aligned_center – array of 2 elements [x,y] that all the model will be centered on for image registration. FIXME: This is the most problematic thing currently, the aligned_center of the model and of the images can be set independently, which will create false results. - In “Load Basis mode”, this parameter is not read, we just use the aligned_center set for the images in the previous fm.klip_dataset and save in basis_filename - In “Save Basis mode”, or “Simple FM mode” we define it and then check that it is the same one used for the images in fm.klip_dataset
mode – deprecated parameter, ignored here and defined in fm.klip_dataset
annuli – deprecated parameter, ignored here and defined in fm.klip_dataset
subsections – deprecated parameter, ignored here and defined in fm.klip_dataset
numthreads – deprecated parameter. All centering are done in fm.klip_dataset
- Returns:
A DiskFM Object
- alloc_fmout(output_img_shape)[source]
Allocates shared memory for the output of the shared memory
- Parameters:
output_img_shape – shape of output image (usually N,y,x,b)
- Returns:
[mp.array to store FM data in, shape of FM data array]
- cleanup_fmout(fmout)[source]
After running KLIP-FM, we need to reshape fmout so that the numKL dimension is the first one and not the last. We also use this function to save the KL basis because it is called by fm.py at the end fm.klip_parallelized
- Parameters:
fmout – numpy array of ouput of FM
- Returns:
Same but cleaned up if necessary
- fm_from_eigen(klmodes=None, evals=None, evecs=None, input_img_shape=None, output_img_shape=None, input_img_num=None, ref_psfs_indicies=None, section_ind=None, aligned_imgs=None, radstart=None, radend=None, phistart=None, phiend=None, padding=None, IOWA=None, ref_center=None, parang=None, numbasis=None, fmout=None, flipx=True, mode=None, **kwargs)[source]
Generate forward models using the KL modes, eigenvectors, and eigenvectors from KLIP. Calls fm.py functions to perform the forward modelling. If we wish to save the KL modes, it save in dictionnaries.
- Parameters:
klmodes – unpertrubed KL modes
evals – eigenvalues of the covariance matrix that generated the KL modes in ascending order(lambda_0 is the 0 index) (shape of [nummaxKL])
evecs – corresponding eigenvectors (shape of [p, nummaxKL])
input_image_shape – 2-D shape of inpt images ([ysize, xsize])
input_img_num – index of sciece frame
ref_psfs_indicies – array of indicies for each reference PSF
section_ind – array indicies into the 2-D x-y image that correspond to this section. Note: needs be called as section_ind[0]
radstart – radius of start of segment
radend – radius of end of segment
phistart – azimuthal start of segment [radians]
phiend – azimuthal end of segment [radians]
padding – amount of padding on each side of sector
IOWA – tuple (IWA,OWA) IWA = Inner working angle & OWA = Outer working angle, both in pixels. It defines the separation interva in which klip will be run.
ref_center – center of image
parang – parallactic angle of input image [DEGREES]
numbasis – array of KL basis cutoffs
fmout – numpy output array for FM output. Shape is (N, y, x, b)
mode – mode of the reduction (‘RDI’, ‘ADI’, ‘SDI’). If RDI only, we only measure the oversubctraction
kwargs – any other variables that we don’t use but are part of the input
- Returns:
None
- fm_parallelized()[source]
Functions like fm.klip_dataset, but it uses previously measured KL modes, section positions, and klip parameter to return the forward modelling. Do not save fits.
- Parameters:
None –
- Returns:
- fmout_np, a numpy array, output of forward modelling
if N_wl = 1, size is [n_KL,x,y]
if N_wl > 1, size is [n_KL,N_wl,x,y]
- load_basis_files(psf_library=None, kl_basis_file=None)[source]
Loads in previously saved basis files and sets variables for fm_from_eigen
- Parameters:
dataset – an instance of Instrument.Data, after fm.klip_dataset. Allow me to pass in the structure some correction parameters set by fm.klip_dataset, such as IWA, OWA, aligned_center. KL basis and sections information are passed via global variables
- Returns:
None
- save_fmout(dataset, fmout, outputdir, fileprefix, numbasis, klipparams=None, calibrate_flux=False, pixel_weights=1, **kwargs)[source]
Uses dataset parameters to save the forward model, the output of fm_paralellized or klip_dataset. No returm, data are saved in “fileprefix” .fits files
- Parameters:
dataset – an instance of Instrument.Data . Will use its dataset.savedata() function to save data
fmout – output of forward modelling.
outputdir – directory to save output files
fileprefix – filename prefix for saved files
numbasis – number of KL basis vectors to use (can be a scalar or list like)
klipparams – string with KLIP-FM parameters
calibrate_flux – if True, flux calibrate the data in the same way as the klipped data
pixel_weights – weights for each pixel for weighted mean. Leave this as a single number for simple mean
- Returns:
None
- save_kl_basis()[source]
Save the KL basis and other needed parameters
- Parameters:
None –
- Returns:
None
- update_disk(model_disk)[source]
Takes model disk and rotates it to the PAs of the input images for use as reference PSFS
The disk can be either an 3D array of shape (wvs,y,x) for data of the same shape or a 2D Array of shape (y,x), in which case, if the dataset is multiwavelength the same model is used for all wavelenths.
- Parameters:
model_disk – Disk to be forward modeled.
- Returns:
None
pyklip.fmlib.extractSpec module
- class pyklip.fmlib.extractSpec.ExtractSpec(inputs_shape, numbasis, sep, pa, input_psfs, input_psfs_wvs, input_psfs_pas=None, datatype='float', stamp_size=None)[source]
Bases:
NoFM
Planet Characterization class. Goal to characterize the astrometry and photometry of a planet
- alloc_fmout(output_img_shape)[source]
Allocates shared memory for the output of the shared memory
- Parameters:
output_img_shape – shape of output image (usually N,y,x,b)
- Returns:
mp.array to store FM data in fmout_shape: shape of FM data array
- Return type:
fmout
- cleanup_fmout(fmout)[source]
After running KLIP-FM, we need to reshape fmout so that the numKL dimension is the first one and not the last
- Parameters:
fmout – numpy array of ouput of FM
- Returns:
same but cleaned up if necessary
- Return type:
fmout
- fm_from_eigen(klmodes=None, evals=None, evecs=None, input_img_shape=None, input_img_num=None, ref_psfs_indicies=None, section_ind=None, section_ind_nopadding=None, aligned_imgs=None, pas=None, wvs=None, radstart=None, radend=None, phistart=None, phiend=None, padding=None, IOWA=None, ref_center=None, parang=None, ref_wv=None, numbasis=None, fmout=None, perturbmag=None, klipped=None, flipx=True, **kwargs)[source]
Generate forward models using the KL modes, eigenvectors, and eigenvectors from KLIP. Calls fm.py functions to perform the forward modelling
- Parameters:
klmodes – unpertrubed KL modes
evals – eigenvalues of the covariance matrix that generated the KL modes in ascending order (lambda_0 is the 0 index) (shape of [nummaxKL])
evecs – corresponding eigenvectors (shape of [p, nummaxKL])
input_image_shape – 2-D shape of inpt images ([ysize, xsize])
input_img_num – index of sciece frame
ref_psfs_indicies – array of indicies for each reference PSF
section_ind – array indicies into the 2-D x-y image that correspond to this section. Note needs be called as section_ind[0]
pas – array of N parallactic angles corresponding to N reference images [degrees]
wvs – array of N wavelengths of those referebce images
radstart – radius of start of segment
radend – radius of end of segment
phistart – azimuthal start of segment [radians]
phiend – azimuthal end of segment [radians]
padding – amount of padding on each side of sector
IOWA – tuple (IWA,OWA) where IWA = Inner working angle and OWA = Outer working angle both in pixels. It defines the separation interva in which klip will be run.
ref_center – center of image
numbasis – array of KL basis cutoffs
parang – parallactic angle of input image [DEGREES]
ref_wv – wavelength of science image
fmout – numpy output array for FM output. Shape is (N, y, x, b)
perturbmag – numpy output for size of linear perturbation. Shape is (N, b)
klipped – PSF subtracted image. Shape of ( size(section), b)
kwargs – any other variables that we don’t use but are part of the input
- generate_models(input_img_shape, section_ind, pas, wvs, radstart, radend, phistart, phiend, padding, ref_center, parang, ref_wv, flipx, stamp_size=None)[source]
Generate model PSFs at the correct location of this segment for each image denoated by its wv and parallactic angle
- Parameters:
pas – array of N parallactic angles corresponding to N images [degrees]
wvs – array of N wavelengths of those images
radstart – radius of start of segment
radend – radius of end of segment
phistart – azimuthal start of segment [radians]
phiend – azimuthal end of segment [radians]
padding – amount of padding on each side of sector
ref_center – center of image
parang – parallactic angle of input image [DEGREES]
ref_wv – wavelength of science image
stamp_size – size of the stamp for spectral extraction
flipx – if True, flip x coordinate in final image
- Returns:
array of size (N, p) where p is the number of pixels in the segment
- Return type:
models
- pyklip.fmlib.extractSpec.invert_spect_fmodel(fmout, dataset, method='JB', units='natural', scaling_factor=1.0)[source]
Greenbaum Nov 2016
- Parameters:
fmout – the forward model matrix which has structure: [numbasis, n_frames, n_frames+1, npix]
dataset – from GPI.GPIData(filelist) – typically set highpass=True also
method – “JB” or “LP” to try the 2 different inversion methods (JB’s or Laurent’s)
units – “natural” means the answer is scaled to the input PSF (default) fmout will be in these units. “scaled” means the output is scaled to “scaling_factor” argument
scaling_factor – multiplies output spectrum and forward model, user set for desired calibration factor. units=”scaled” must be set in args for this to work!
- Returns:
A tuple containing the spectrum and the forward model (spectrum, forwardmodel) spectrum shape:(len(numbasis), nwav)
pyklip.fmlib.fmpsf module
- class pyklip.fmlib.fmpsf.FMPlanetPSF(inputs_shape, numbasis, sep, pa, dflux, input_psfs, input_wvs, flux_conversion=None, spectrallib=None, spectrallib_units='flux', star_spt=None, refine_fit=False, field_dependent_correction=None, input_psfs_pas=None)[source]
Bases:
NoFM
Forward models the PSF of the planet through KLIP. Returns the forward modelled planet PSF
- alloc_fmout(output_img_shape)[source]
Allocates shared memory for the output of the shared memory
- Parameters:
output_img_shape – shape of output image (usually N,y,x,b)
- Returns:
mp.array to store FM data in fmout_shape: shape of FM data array
- Return type:
fmout
- alloc_perturbmag(output_img_shape, numbasis)[source]
Allocates shared memory to store the fractional magnitude of the linear KLIP perturbation Stores a number for each frame = max(oversub + selfsub)/std(PCA(image))
- Parameters:
output_img_shape – shape of output image (usually N,y,x,b)
numbasis – array/list of number of KL basis cutoffs requested
- Returns:
mp.array to store linaer perturbation magnitude perturbmag_shape: shape of linear perturbation magnitude
- Return type:
perturbmag
- cleanup_fmout(fmout)[source]
After running KLIP-FM, we need to reshape fmout so that the numKL dimension is the first one and not the last
- Parameters:
fmout – numpy array of ouput of FM
- Returns:
same but cleaned up if necessary
- Return type:
fmout
- fm_from_eigen(klmodes=None, evals=None, evecs=None, input_img_shape=None, input_img_num=None, ref_psfs_indicies=None, section_ind=None, aligned_imgs=None, pas=None, wvs=None, radstart=None, radend=None, phistart=None, phiend=None, padding=None, IOWA=None, ref_center=None, parang=None, ref_wv=None, numbasis=None, fmout=None, perturbmag=None, klipped=None, covar_files=None, flipx=True, rdi_psfs=None, **kwargs)[source]
Generate forward models using the KL modes, eigenvectors, and eigenvectors from KLIP. Calls fm.py functions to perform the forward modelling
- Parameters:
klmodes – unpertrubed KL modes
evals – eigenvalues of the covariance matrix that generated the KL modes in ascending order (lambda_0 is the 0 index) (shape of [nummaxKL])
evecs – corresponding eigenvectors (shape of [p, nummaxKL])
input_image_shape – 2-D shape of inpt images ([ysize, xsize])
input_img_num – index of sciece frame
ref_psfs_indicies – array of indicies for each reference PSF
section_ind – array indicies into the 2-D x-y image that correspond to this section. Note needs be called as section_ind[0]
pas – array of N parallactic angles corresponding to N reference images [degrees]
wvs – array of N wavelengths of those referebce images
radstart – radius of start of segment
radend – radius of end of segment
phistart – azimuthal start of segment [radians]
phiend – azimuthal end of segment [radians]
padding – amount of padding on each side of sector
IOWA – tuple (IWA,OWA) where IWA = Inner working angle and OWA = Outer working angle both in pixels. It defines the separation interva in which klip will be run.
ref_center – center of image
numbasis – array of KL basis cutoffs
parang – parallactic angle of input image [DEGREES]
ref_wv – wavelength of science image
fmout – numpy output array for FM output. Shape is (N, y, x, b)
perturbmag – numpy output for size of linear perturbation. Shape is (N, b)
klipped – PSF subtracted image. Shape of ( size(section), b)
ref_rdi_indices – array of N+M indicating N ADI/SDI images (0’s) and M RDI images (1;s).
rdi_psfs – array of M RDI reference images in this sector.
kwargs – any other variables that we don’t use but are part of the input
- generate_models(input_img_shape, section_ind, pas, wvs, radstart, radend, phistart, phiend, padding, ref_center, parang, ref_wv, flipx, rdi_indices)[source]
Generate model PSFs at the correct location of this segment for each image denoated by its wv and parallactic angle
- Parameters:
pas – array of N parallactic angles corresponding to N images [degrees]
wvs – array of N wavelengths of those images
radstart – radius of start of segment
radend – radius of end of segment
phistart – azimuthal start of segment [radians]
phiend – azimuthal end of segment [radians]
padding – amount of padding on each side of sector
ref_center – center of image
parang – parallactic angle of input image [DEGREES]
ref_wv – wavelength of science image
flipx – if True, flip x coordinate in final image
rdi_indices – array of N corresponding to whether it is (1) or isn’t (0) an RDI frame
- Returns:
array of size (N, p) where p is the number of pixels in the segment
- Return type:
models
- save_fmout(dataset, fmout, outputdir, fileprefix, numbasis, klipparams=None, calibrate_flux=False, spectrum=None, pixel_weights=1)[source]
Saves the FM planet PSFs to disk. Saves both a KL mode cube and spectral cubes
- Parameters:
dataset – Instruments.Data instance. Will use its dataset.savedata() function to save data
fmout – the fmout data passed from fm.klip_parallelized which is passed as the output of cleanup_fmout
outputdir – output directory
fileprefix – the fileprefix to prepend the file name
numbasis – KL mode cutoffs used
klipparams – string with KLIP-FM parameters
calibrate_flux – if True, flux calibrate the data in the same way as the klipped data
spectrum – if not None, spectrum to weight the data by. Length same as dataset.wvs
pixel_weights – weights for each pixel for weighted mean. Leave this as a single number for simple mean
pyklip.fmlib.matchedFilter module
- class pyklip.fmlib.matchedFilter.MatchedFilter(inputs_shape, numbasis, input_psfs, input_psfs_wvs, spectrallib, save_per_sector=None, datatype='float', fakes_sepPa_list=None, disable_FM=None, true_fakes_pos=None, ref_center=None, flipx=None, rm_edge=None, edge_threshold=None, planet_radius=None, background_width=None, save_bbfm=None, save_fm=None, save_fmout_path=None)[source]
Bases:
NoFM
Matched filter with forward modelling.
- alloc_fmout(output_img_shape)[source]
Allocates shared memory for the output of the shared memory
- Parameters:
output_img_shape – Not used
- Returns:
mp.array to store auxilliary data in fmout_shape: shape of auxilliary array = (3,self.N_spectra,self.N_numbasis,self.N_frames,self.ny,self.nx)
- The 3 is for saving the different term of the matched filter:
0: dot product 1: square of the norm of the model 2: Local estimated variance of the data 3: Number of pixels used in the matched filter
- Return type:
fmout
- fm_end_sector(interm_data=None, fmout=None, sector_index=None, section_indicies=None)[source]
Save the fmout object at the end of each sector if save_per_sector was defined when initializing the class.
- fm_from_eigen(klmodes=None, evals=None, evecs=None, input_img_shape=None, input_img_num=None, ref_psfs_indicies=None, section_ind=None, section_ind_nopadding=None, aligned_imgs=None, pas=None, wvs=None, radstart=None, radend=None, phistart=None, phiend=None, padding=None, IOWA=None, ref_center=None, parang=None, ref_wv=None, numbasis=None, fmout=None, perturbmag=None, klipped=None, flipx=True, rdi_psfs=None, **kwargs)[source]
Calculate and project the FM at every pixel of the sector. Store the result in fmout.
- Parameters:
klmodes – unpertrubed KL modes
evals – eigenvalues of the covariance matrix that generated the KL modes in ascending order (lambda_0 is the 0 index) (shape of [nummaxKL])
evecs – corresponding eigenvectors (shape of [p, nummaxKL])
input_img_shape – 2-D shape of inpt images ([ysize, xsize])
input_img_num – index of sciece frame
ref_psfs_indicies – array of indicies for each reference PSF
section_ind – array indicies into the 2-D x-y image that correspond to this section. Note needs be called as section_ind[0]
pas – array of N parallactic angles corresponding to N reference images [degrees]
wvs – array of N wavelengths of those referebce images
radstart – radius of start of segment
radend – radius of end of segment
phistart – azimuthal start of segment [radians]
phiend – azimuthal end of segment [radians]
padding – amount of padding on each side of sector
IOWA – tuple (IWA,OWA) where IWA = Inner working angle and OWA = Outer working angle both in pixels. It defines the separation interva in which klip will be run.
ref_center – center of image
numbasis – array of KL basis cutoffs
parang – parallactic angle of input image [DEGREES]
ref_wv – wavelength of science image
fmout – numpy output array for FM output. Shape is (N, y, x, b)
klipped – array of shape (p,b) that is the PSF subtracted data for each of the b KLIP basis cutoffs. If numbasis was an int, then sub_img_row_selected is just an array of length p
rdi_psfs – array of M RDI reference images in this sector.
kwargs – any other variables that we don’t use but are part of the input
- generate_model_sci(input_img_shape, section_ind, pa, wv, radstart, radend, phistart, phiend, padding, ref_center, parang, ref_wv, sep_fk, pa_fk, flipx)[source]
Generate model PSFs at the correct location of this segment of the science image denotated by its wv and parallactic angle.
- Parameters:
input_img_shape – 2-D shape of inpt images ([ysize, xsize])
section_ind – array indicies into the 2-D x-y image that correspond to this section. Note needs be called as section_ind[0]
pa – parallactic angle of the science image [degrees]
wv – wavelength of the science image
radstart – radius of start of segment (not used)
radend – radius of end of segment (not used)
phistart – azimuthal start of segment [radians] (not used)
phiend – azimuthal end of segment [radians] (not used)
padding – amount of padding on each side of sector
ref_center – center of image
parang – parallactic angle of input image [DEGREES] (not used)
ref_wv – wavelength of science image
sep_fk – separation of the planet to be injected.
pa_fk – position angle of the planet to be injected.
flipx – if True, flip x coordinate in final image
- Return: (models, mask)
models: vector of size p where p is the number of pixels in the segment mask: vector of size p where p is the number of pixels in the segment
if pixel == 1: arc shape where to calculate the standard deviation if pixel == 2: 7 pixels disk around the position of the planet.
- generate_models(input_img_shape, section_ind, pas, wvs, radstart, radend, phistart, phiend, padding, ref_center, parang, ref_wv, sep_fk, pa_fk, flipx, rdi_indices)[source]
Generate model PSFs at the correct location of this segment for each image denotated by its wv and parallactic angle.
- Parameters:
input_img_shape – 2-D shape of inpt images ([ysize, xsize])
section_ind – array indicies into the 2-D x-y image that correspond to this section. Note needs be called as section_ind[0]
pas – array of N parallactic angles corresponding to N images [degrees]
wvs – array of N wavelengths of those images
radstart – radius of start of segment (not used)
radend – radius of end of segment (not used)
phistart – azimuthal start of segment [radians] (not used)
phiend – azimuthal end of segment [radians] (not used)
padding – amount of padding on each side of sector
ref_center – center of image
parang – parallactic angle of input image [DEGREES] (not used)
ref_wv – wavelength of science image
sep_fk – separation of the planet to be injected.
pa_fk – position angle of the planet to be injected.
flipx – if True, flip x coordinate in final image
- Returns:
array of size (N, p) where p is the number of pixels in the segment
- Return type:
models
- save_fmout(dataset, fmout, outputdir, fileprefix, numbasis, klipparams=None, calibrate_flux=False, spectrum=None, pixel_weights=1)[source]
Saves the fmout data to disk following the instrument’s savedata function
- Parameters:
dataset – Instruments.Data instance. Will use its dataset.savedata() function to save data
fmout – the fmout data passed from fm.klip_parallelized which is passed as the output of cleanup_fmout
outputdir – output directory
fileprefix – the fileprefix to prepend the file name
numbasis – KL mode cutoffs used
klipparams – string with KLIP-FM parameters
calibrate_flux – if True, flux calibrate the data (if applicable)
spectrum – if not None, the spectrum to weight the data by. Length same as dataset.wvs
pixel_weights – weights for each pixel for weighted mean. [Not used in matched filter]
- skip_section(radstart, radend, phistart, phiend, flipx=True)[source]
Returns a boolean indicating if the section defined by (radstart, radend, phistart, phiend) should be skipped. When True is returned the current section in the loop in klip_parallelized() is skipped.
- Parameters:
radstart – minimum radial distance of sector [pixels]
radend – maximum radial distance of sector [pixels]
phistart – minimum azimuthal coordinate of sector [radians]
phiend – maximum azimuthal coordinate of sector [radians]
flipx – if True, flip x coordinate in final image
- Returns:
False so by default it never skips.
- Return type:
Boolean
- pyklip.fmlib.matchedFilter.calculate_fm_opti(delta_KL, original_KL, sci, model_sci_fk, delta_KL_fk, original_KL_fk)[source]
Optimized version for calculate_fm() (if numbas) for a single numbasis.
Calculate what the PSF looks up post-KLIP using knowledge of the input PSF, assumed spectrum of the science target, and the partially calculated KL modes (Delta Z_k^lambda in Laurent’s paper). If inputflux is None, the spectral dependence has already been folded into delta_KL_nospec (treat it as delta_KL).
Note: if inputflux is None and delta_KL_nospec has three dimensions (ie delta_KL_nospec was calculated using pertrurb_nospec() or perturb_nospec_modelsBased()) then only klipped_oversub and klipped_selfsub are returned. Besides they will have an extra first spectral dimension.
- Parameters:
delta_KL_nospec – perturbed KL modes but without the spectral info. delta_KL = spectrum x delta_Kl_nospec. Shape is (numKL, wv, pix). If inputflux is None, delta_KL_nospec = delta_KL
orignal_KL – unpertrubed KL modes (array of size [numbasis, numpix])
sci – array of size p representing the science data
model_sci – array of size p corresponding to the PSF of the science frame
input_spectrum – array of size wv with the assumed spectrum of the model
If delta_KL_nospec does NOT include a spectral dimension or if inputflux is not None: :returns:
- array of shape (b,p) showing the forward modelled PSF
Skipped if inputflux = None, and delta_KL_nospec has 3 dimensions.
klipped_oversub: array of shape (b, p) showing the effect of oversubtraction as a function of KL modes klipped_selfsub: array of shape (b, p) showing the effect of selfsubtraction as a function of KL modes Note: psf_FM = model_sci - klipped_oversub - klipped_selfsub to get the FM psf as a function of K Lmodes
(shape of b,p)
- Return type:
fm_psf
If inputflux = None and if delta_KL_nospec include a spectral dimension: :returns: Sum(<S|KL>KL) with klipped_oversub.shape = (1,Npix)
klipped_selfsub: Sum(<N|DKL>KL) + Sum(<N|KL>DKL) with klipped_selfsub.shape = (1,N_lambda or N_ref,N_pix)
- Return type:
klipped_oversub
pyklip.fmlib.nofm module
- class pyklip.fmlib.nofm.NoFM(inputs_shape, numbasis)[source]
Bases:
object
Super class for all forward modelling classes. Has fall-back functions for all fm dependent calls so that each FM class does not need to implement functions it doesn’t want to. Should do no forward modelling and just do regular KLIP by itself
- alloc_fmout(output_img_shape)[source]
Allocates shared memory for the output of the forward modelling
- Parameters:
output_img_shape – shape of output image (usually N,y,x,b)
- Returns:
mp.array to store FM data in fmout_shape: shape of FM data array
- Return type:
fmout
- alloc_interm(max_sector_size, numsciframes)[source]
Allocates shared memory array for intermediate step
Intermediate step is allocated for a sector by sector basis
- Parameters:
max_sector_size – number of pixels in this sector. Max because this can be variable. Stupid rotating sectors
- Returns:
mp.array to store intermediate products from one sector in interm_shape:shape of interm array (used to convert to numpy arrays)
- Return type:
interm
- alloc_output()[source]
Allocates shared memory array for final output
Only use multiprocessing data structors as we are using the multiprocessing class
Args:
- Returns:
mp.array to store final outputs in (shape of (N*wv, y, x, numbasis)) outputs_shape: shape of outputs array to convert to numpy arrays
- Return type:
outputs
- alloc_perturbmag(output_img_shape, numbasis)[source]
Allocates shared memory to store the fractional magnitude of the linear KLIP perturbation
- Parameters:
output_img_shape – shape of output image (usually N,y,x,b)
numbasis – array/list of number of KL basis cutoffs requested
- Returns:
mp.array to store linaer perturbation magnitude perturbmag_shape: shape of linear perturbation magnitude
- Return type:
perturbmag
- cleanup_fmout(fmout)[source]
After running KLIP-FM, if there’s anything to do to the fmout array (such as reshaping it), now’s the time to do that before outputting it
- Parameters:
fmout – numpy array of ouput of FM
- Returns:
same but cleaned up if necessary
- Return type:
fmout
- fm_end_sector(**kwargs)[source]
Does some forward modelling at the end of a sector after all images have been klipped for that sector.
- fm_from_eigen(**kwargs)[source]
Generate forward models using the KL modes, eigenvectors, and eigenvectors from KLIP This is called immediately after regular KLIP
- save_fmout(dataset, fmout, outputdir, fileprefix, numbasis, klipparams=None, calibrate_flux=False, spectrum=None, pixel_weights=1)[source]
Saves the fmout data to disk following the instrument’s savedata function
- Parameters:
dataset – Instruments.Data instance. Will use its dataset.savedata() function to save data
fmout – the fmout data passed from fm.klip_parallelized which is passed as the output of cleanup_fmout
outputdir – output directory
fileprefix – the fileprefix to prepend the file name
numbasis – KL mode cutoffs used
klipparams – string with KLIP-FM parameters
calibrate_flux – if True, flux calibrate the data (if applicable)
spectrum – if not None, the spectrum to weight the data by. Length same as dataset.wvs
pixel_weights – weights for each pixel for weighted mean. Leave this as a single number for simple mean
- skip_section(radstart, radend, phistart, phiend, flipx=None)[source]
Returns a boolean indicating if the section defined by (radstart, radend, phistart, phiend) should be skipped. When True is returned the current section in the loop in klip_parallelized() is skipped.
- Parameters:
radstart – minimum radial distance of sector [pixels]
radend – maximum radial distance of sector [pixels]
phistart – minimum azimuthal coordinate of sector [radians]
phiend – maximum azimuthal coordinate of sector [radians]
flipx – if True, flip x coordinate in final image
- Returns:
False so by default it never skips.
- Return type:
Boolean