lenstronomy.Workflow package

Submodules

lenstronomy.Workflow.alignment_matching module

class lenstronomy.Workflow.alignment_matching.AlignmentFitting(multi_band_list, kwargs_model, kwargs_params, band_index=0, likelihood_mask_list=None)[source]

Bases: object

class which executes the different sampling methods

pso(n_particles=10, n_iterations=10, lowerLimit=-0.2, upperLimit=0.2, threadCount=1, mpi=False, print_key='default')[source]

returns the best fit for the lense model on catalogue basis with particle swarm optimizer

class lenstronomy.Workflow.alignment_matching.AlignmentLikelihood(multi_band_list, kwargs_model, kwargs_params, band_index=0, likelihood_mask_list=None)[source]

Bases: object

computeLikelihood(ctx)[source]
static get_args(kwargs_data)[source]
Parameters:kwargs_data
Returns:
likelihood(a)[source]
num_param
setup()[source]
update_data(args)[source]
Parameters:args
Returns:
update_multi_band(args)[source]
Parameters:args – list of parameters
Returns:updated multi_band_list

lenstronomy.Workflow.fitting_sequence module

class lenstronomy.Workflow.fitting_sequence.FittingSequence(kwargs_data_joint, kwargs_model, kwargs_constraints, kwargs_likelihood, kwargs_params, mpi=False, verbose=True)[source]

Bases: object

class to define a sequence of fitting applied, inherit the Fitting class this is a Workflow manager that allows to update model configurations before executing another step in the modelling The user can take this module as an example of how to create their own workflows or build their own around the FittingSequence

align_images(n_particles=10, n_iterations=10, lowerLimit=-0.2, upperLimit=0.2, threadCount=1, compute_bands=None)[source]

aligns the coordinate systems of different exposures within a fixed model parameterisation by executing a PSO with relative coordinate shifts as free parameters

Parameters:
  • n_particles – number of particles in the Particle Swarm Optimization
  • n_iterations – number of iterations in the optimization process
  • lowerLimit – lower limit of relative shift
  • upperLimit – upper limit of relative shift
  • verbose – bool, print statements
  • compute_bands – bool list, if multiple bands, this process can be limited to a subset of bands
Returns:

best_fit(bijective=False)[source]
Parameters:bijective – bool, if True, the mapping of image2source_plane and the mass_scaling parameterisation are inverted. If you do not use those options, there is no effect.
Returns:best fit model of the current state of the FittingSequence class
best_fit_likelihood

returns the log likelihood of the best fit model of the current state of this class

Returns:log likelihood, float
bic

returns the bayesian information criterion of the model. :return: bic value, float

fit_sequence(fitting_list)[source]
Parameters:fitting_list – list of [[‘string’, {kwargs}], ..] with ‘string being the specific fitting option and kwargs being the arguments passed to this option
Returns:fitting results
fix_not_computed(free_bands)[source]

fixes lens model parameters of imaging bands/frames that are not computed and frees the parameters of the other lens models to the initial kwargs_fixed options

Parameters:free_bands – bool list of length of imaging bands in order of imaging bands, if False: set fixed lens model
Returns:None
kwargs_fixed()[source]

returns the updated kwargs_fixed from the update Manager

Returns:list of fixed kwargs, see UpdateManager()
likelihoodModule
Returns:Likelihood() class instance reflecting the current state of Fittingsequence
mcmc(n_burn, n_run, walkerRatio, sigma_scale=1, threadCount=1, init_samples=None, re_use_samples=True, sampler_type='EMCEE', progress=True, backup_filename=None, start_from_backup=False)[source]

MCMC routine

Parameters:
  • n_burn – number of burn in iterations (will not be saved)
  • n_run – number of MCMC iterations that are saved
  • walkerRatio – ratio of walkers/number of free parameters
  • sigma_scale – scaling of the initial parameter spread relative to the width in the initial settings
  • threadCount – number of CPU threads. If MPI option is set, threadCount=1
  • init_samples – initial sample from where to start the MCMC process
  • re_use_samples – bool, if True, re-uses the samples described in init_samples.nOtherwise starts from scratch.
  • sampler_type – string, which MCMC sampler to be used. Options are: ‘EMCEE’
  • progress – boolean, if True shows progress bar in EMCEE
Returns:

list of output arguments, e.g. MCMC samples, parameter names, logL distances of all samples specified by the specific sampler used

nested_sampling(sampler_type='MULTINEST', kwargs_run={}, prior_type='uniform', width_scale=1, sigma_scale=1, output_basename='chain', remove_output_dir=True, dypolychord_dynamic_goal=0.8, polychord_settings={}, dypolychord_seed_increment=200, output_dir='nested_sampling_chains', dynesty_bound='multi', dynesty_sample='auto')[source]

Run (Dynamic) Nested Sampling algorithms, depending on the type of algorithm.

Parameters:
  • sampler_type – ‘MULTINEST’, ‘DYPOLYCHORD’, ‘DYNESTY’
  • kwargs_run – keywords passed to the core sampling method
  • prior_type – ‘uniform’ of ‘gaussian’, for converting the unit hypercube to param cube
  • width_scale – scale the width (lower/upper limits) of the parameters space by this factor
  • sigma_scale – if prior_type is ‘gaussian’, scale the gaussian sigma by this factor
  • output_basename – name of the folder in which the core MultiNest/PolyChord code will save output files
  • remove_output_dir – if True, the above folder is removed after completion
  • dypolychord_dynamic_goal – dynamic goal for DyPolyChord (trade-off between evidence (0) and posterior (1) computation)
  • polychord_settings – settings dictionary to send to pypolychord. Check dypolychord documentation for details.
  • dypolychord_seed_increment – seed increment for dypolychord with MPI. Check dypolychord documentation for details.
  • dynesty_bound – see https://dynesty.readthedocs.io for details
  • dynesty_sample – see https://dynesty.readthedocs.io for details
Returns:

list of output arguments : samples, mean inferred values, log-likelihood, log-evidence, error on log-evidence for each sample

param_class
Returns:Param() class instance reflecting the current state of Fittingsequence
psf_iteration(compute_bands=None, **kwargs_psf_iter)[source]

iterative PSF reconstruction

Parameters:
  • compute_bands – bool list, if multiple bands, this process can be limited to a subset of bands
  • kwargs_psf_iter – keyword arguments as used or available in PSFIteration.update_iterative() definition
Returns:

0, updated PSF is stored in self.multi_band_list

pso(n_particles, n_iterations, sigma_scale=1, print_key='PSO', threadCount=1)[source]

Particle Swarm Optimization

Parameters:
  • n_particles – number of particles in the Particle Swarm Optimization
  • n_iterations – number of iterations in the optimization process
  • sigma_scale – scaling of the initial parameter spread relative to the width in the initial settings
  • print_key – string, printed text when executing this routine
  • threadCount – number of CPU threads. If MPI option is set, threadCount=1
Returns:

result of the best fit, the chain of the best fit parameter after each iteration, list of parameters in same order

set_param_value(**kwargs)[source]

Set a parameter to a specific value. kwargs are below. :param lens: [[i_model, [‘param1’, ‘param2’,…], […]] :type lens: :param source: [[i_model, [‘param1’, ‘param2’,…], […]] :type source: :param lens_light: [[i_model, [‘param1’, ‘param2’,…], […]] :type lens_light: :param ps: [[i_model, [‘param1’, ‘param2’,…], […]] :type ps: :return: 0, the value of the param is overwritten :rtype:

simplex(n_iterations, method='Nelder-Mead')[source]

Downhill simplex optimization using the Nelder-Mead algorithm.

Parameters:
  • n_iterations – maximum number of iterations to perform
  • method – the optimization method used, see documentation in scipy.optimize.minimize
Returns:

result of the best fit

update_settings(kwargs_model={}, kwargs_constraints={}, kwargs_likelihood={}, lens_add_fixed=[], source_add_fixed=[], lens_light_add_fixed=[], ps_add_fixed=[], cosmo_add_fixed=[], lens_remove_fixed=[], source_remove_fixed=[], lens_light_remove_fixed=[], ps_remove_fixed=[], cosmo_remove_fixed=[], change_source_lower_limit=None, change_source_upper_limit=None)[source]

updates lenstronomy settings “on the fly”

Parameters:
  • kwargs_model – kwargs, specified keyword arguments overwrite the existing ones
  • kwargs_constraints – kwargs, specified keyword arguments overwrite the existing ones
  • kwargs_likelihood – kwargs, specified keyword arguments overwrite the existing ones
  • lens_add_fixed – [[i_model, [‘param1’, ‘param2’,…], […]]
  • source_add_fixed – [[i_model, [‘param1’, ‘param2’,…], […]]
  • lens_light_add_fixed – [[i_model, [‘param1’, ‘param2’,…], […]]
  • ps_add_fixed – [[i_model, [‘param1’, ‘param2’,…], […]]
  • cosmo_add_fixed – [‘param1’, ‘param2’,…]
  • lens_remove_fixed – [[i_model, [‘param1’, ‘param2’,…], […]]
  • source_remove_fixed – [[i_model, [‘param1’, ‘param2’,…], […]]
  • lens_light_remove_fixed – [[i_model, [‘param1’, ‘param2’,…], […]]
  • ps_remove_fixed – [[i_model, [‘param1’, ‘param2’,…], […]]
  • cosmo_remove_fixed – [‘param1’, ‘param2’,…]
Returns:

0, the settings are overwritten for the next fitting step to come

update_state(kwargs_update)[source]

updates current best fit state to the input model keywords specified.

Parameters:kwargs_update – format of kwargs_result
Returns:None

lenstronomy.Workflow.psf_fitting module

class lenstronomy.Workflow.psf_fitting.PsfFitting(image_model_class)[source]

Bases: object

class to find subsequently a better psf The method make use of a model and subtracts all the non-point source components of the model from the data. If the model is sufficient, then the data will be a (better) representation of the actual PSF. The method cuts out those point sources and combines them to update the estimate of the PSF. This is done in an iterative procedure as the model components of the extended features is PSF-dependent (hopefully not too much).

Various options can be chosen. There is no guarantee that the method works for specific data and models.

‘stacking_method’: ‘median’, ‘mean’; the different estimates of the PSF are stacked and combined together. The choices are:
‘mean’: mean of pixel values as the estimator (not robust to outliers) ‘median’: median of pixel values as the estimator (outlier rejection robust but needs >2 point sources in the data
‘block_center_neighbour’: angle, radius of neighbouring point sources around their centers the estimates is ignored.
Default is zero, meaning a not optimal subtraction of the neighbouring point sources might contaminate the estimate.
‘keep_error_map’: bool, if True, does not replace the error term associated with the PSF estimate.
If false, re-estimates the variance between the PSF estimates.
‘psf_symmetry’: number of rotational invariant symmetries in the estimated PSF.
=1 mean no additional symmetries. =4 means 90 deg symmetry. This is enforced by a rotatioanl stack according to the symmetry specified. These additional imposed symmetries can help stabelize the PSF estimate when there are limited constraints/number of point sources in the image.

The procedure only requires and changes the ‘point_source_kernel’ in the PSF() class and the ‘psf_error_map’. Any previously set subgrid kernels or pixel_kernels are removed and constructed from the ‘point_source_kernel’.

static combine_psf(kernel_list_new, kernel_old, factor=1.0, stacking_option='median', symmetry=1)[source]

updates psf estimate based on old kernel and several new estimates :param kernel_list_new: list of new PSF kernels estimated from the point sources in the image (un-normalized) :param kernel_old: old PSF kernel :param factor: weight of updated estimate based on new and old estimate, factor=1 means new estimate, factor=0 means old estimate :param stacking_option: option of stacking, mean or median :param symmetry: imposed symmetry of PSF estimate :return: updated PSF estimate and error_map associated with it

cutout_psf(ra_image, dec_image, x, y, image_list, kernelsize, kernel_init, block_center_neighbour=0)[source]
Parameters:
  • x
  • y
  • image_list – list of images (i.e. data - all models subtracted, except a single point source)
  • kernelsize
Returns:

static cutout_psf_single(x, y, image, mask, kernelsize, kernel_init)[source]
Parameters:
  • x – x-coordinate of point soure
  • y – y-coordinate of point source
  • image – image (i.e. data - all models subtracted, except a single point source)
  • mask – mask of pixels in the image not to be considered in the PSF estimate (being replaced by kernel_init)
  • kernelsize – width in pixel of the kernel
  • kernel_init – initial guess of kernel (pixels that are masked are replaced by those values)
Returns:

estimate of the PSF based on the image and position of the point source

error_map_estimate(kernel, star_cutout_list, amp, x_pos, y_pos, error_map_radius=None, block_center_neighbour=0)[source]

provides a psf_error_map based on the goodness of fit of the given PSF kernel on the point source cutouts, their estimated amplitudes and positions

Parameters:
  • kernel – PSF kernel
  • star_cutout_list – list of 2d arrays of cutouts of the point sources with all other model components subtracted
  • amp – list of amplitudes of the estimated PSF kernel
  • x_pos – pixel position (in original data unit, not in cutout) of the point sources (same order as amp and star cutouts)
  • y_pos – pixel position (in original data unit, not in cutout) of the point sources (same order as amp and star cutouts)
  • error_map_radius – float, radius (in arc seconds) of the outermost error in the PSF estimate (e.g. to avoid double counting of overlapping PSF erros)
  • block_center_neighbour – angle, radius of neighbouring point sources around their centers the estimates is ignored. Default is zero, meaning a not optimal subtraction of the neighbouring point sources might contaminate the estimate.
Returns:

relative uncertainty in the psf model (in quadrature) per pixel based on residuals achieved in the image

image_single_point_source(image_model_class, kwargs_params)[source]

return model without including the point source contributions as a list (for each point source individually)

Parameters:
  • image_model_class – ImageModel class instance
  • kwargs_params – keyword arguments of model component keyword argument lists
Returns:

list of images with point source isolated

static mask_point_source(x_pos, y_pos, x_grid, y_grid, radius, i=0)[source]
Parameters:
  • x_pos – x-position of list of point sources
  • y_pos – y-position of list of point sources
  • x_grid – x-coordinates of grid
  • y_grid – y-coordinates of grid
  • i – index of point source not to mask out
  • radius – radius to mask out other point sources
Returns:

a mask of the size of the image with cutouts around the position

update_iterative(kwargs_psf, kwargs_params, num_iter=10, keep_psf_error_map=True, no_break=True, verbose=True, **kwargs_psf_update)[source]
Parameters:
  • kwargs_psf – keyword arguments to construct the PSF() class
  • kwargs_params – keyword arguments of the parameters of the model components (e.g. ‘kwargs_lens’ etc)
  • num_iter – number of iterations in the PSF fitting and image fitting process
  • keep_psf_error_map – boolean, if True keeps previous psf_error_map
  • no_break – boolean, if True, runs until the end regardless of the next step getting worse, and then reads out the overall best fit
  • verbose – print statements informing about progress of iterative procedure
  • kwargs_psf_update – keyword arguments providing the settings for a single iteration of the PSF, as being passed to update_psf() method
Returns:

keyword argument of PSF constructor for PSF() class with updated PSF

update_psf(kwargs_psf, kwargs_params, stacking_method='median', psf_symmetry=1, psf_iter_factor=0.2, block_center_neighbour=0, error_map_radius=None, block_center_neighbour_error_map=None)[source]
Parameters:
  • kwargs_psf – keyword arguments to construct the PSF() class
  • kwargs_params – keyword arguments of the parameters of the model components (e.g. ‘kwargs_lens’ etc)
  • stacking_method – ‘median’, ‘mean’; the different estimates of the PSF are stacked and combined together. The choices are: ‘mean’: mean of pixel values as the estimator (not robust to outliers) ‘median’: median of pixel values as the estimator (outlier rejection robust but needs >2 point sources in the data
  • psf_symmetry – number of rotational invariant symmetries in the estimated PSF. =1 mean no additional symmetries. =4 means 90 deg symmetry. This is enforced by a rotatioanl stack according to the symmetry specified. These additional imposed symmetries can help stabelize the PSF estimate when there are limited constraints/number of point sources in the image.
  • psf_iter_factor – factor in (0, 1] of ratio of old vs new PSF in the update in the iteration.
  • block_center_neighbour – angle, radius of neighbouring point sources around their centers the estimates is ignored. Default is zero, meaning a not optimal subtraction of the neighbouring point sources might contaminate the estimate.
  • block_center_neighbour_error_map – angle, radius of neighbouring point sources around their centers the estimates of the ERROR MAP is ignored. If None, then the value of block_center_neighbour is used (recommended)
  • error_map_radius – float, radius (in arc seconds) of the outermost error in the PSF estimate (e.g. to avoid double counting of overlapping PSF errors), if None, all of the pixels are considered (unless blocked through other means)
Returns:

kwargs_psf_new, logL_after, error_map

lenstronomy.Workflow.update_manager module

class lenstronomy.Workflow.update_manager.UpdateManager(kwargs_model, kwargs_constraints, kwargs_likelihood, kwargs_params)[source]

Bases: object

this class manages the parameter constraints as they may evolve through the steps of the modeling. This includes: keeping certain parameters fixed during one modelling step

best_fit(bijective=False)[source]
fix_image_parameters(image_index=0)[source]

fixes all parameters that are only assigned to a specific image. This allows to sample only parameters that constraint by the fitting of a sub-set of the images.

Parameters:image_index – index
Returns:
fixed_kwargs
init_kwargs
param_class
Returns:instance of the Param class with the recent options and bounds
parameter_state
Returns:parameter state saved in this class
set_init_state()[source]

set the current state of the parameters to the initial one.

Returns:
sigma_kwargs
update_fixed(lens_add_fixed=[], source_add_fixed=[], lens_light_add_fixed=[], ps_add_fixed=[], special_add_fixed=[], lens_remove_fixed=[], source_remove_fixed=[], lens_light_remove_fixed=[], ps_remove_fixed=[], special_remove_fixed=[])[source]

adds the values of the keyword arguments that are stated in the _add_fixed to the existing fixed arguments.

Parameters:
  • lens_add_fixed
  • source_add_fixed
  • lens_light_add_fixed
  • ps_add_fixed
  • special_add_fixed
Returns:

updated kwargs fixed

update_limits(change_source_lower_limit=None, change_source_upper_limit=None)[source]

updates the limits (lower and upper) of the update manager instance

Parameters:change_source_lower_limit – [[i_model, [‘param_name’, …], [value1, value2, …]]]
Returns:updates internal state of lower and upper limits accessible from outside
update_options(kwargs_model, kwargs_constraints, kwargs_likelihood)[source]

updates the options by overwriting the kwargs with the new ones being added/changed WARNING: some updates may not be valid depending on the model options. Use carefully!

Parameters:
  • kwargs_model
  • kwargs_constraints
  • kwargs_likelihood
Returns:

update_param_state(kwargs_lens=None, kwargs_source=None, kwargs_lens_light=None, kwargs_ps=None, kwargs_special=None, kwargs_extinction=None)[source]

updates the temporary state of the parameters being saved. ATTENTION: Any previous knowledge gets lost if you call this function

Parameters:
  • kwargs_lens
  • kwargs_source
  • kwargs_lens_light
  • kwargs_ps
  • kwargs_special
  • kwargs_extinction
Returns:

update_param_value(lens=[], source=[], lens_light=[], ps=[])[source]

Set a model parameter to a specific value. :param lens: [[i_model, [‘param1’, ‘param2’,…], […]] :type lens: :param source: [[i_model, [‘param1’, ‘param2’,…], […]] :type source: :param lens_light: [[i_model, [‘param1’, ‘param2’,…], […]] :type lens_light: :param ps: [[i_model, [‘param1’, ‘param2’,…], […]] :type ps: :return: 0, the value of the param is overwritten :rtype:

Module contents