make_model_image#
- photutils.datasets.make_model_image(shape, model, params_table, *, model_shape=None, bbox_factor=None, x_name='x_0', y_name='y_0', discretize_method='center', discretize_oversample=10, progress_bar=False)[source]#
Make a 2D image containing sources generated from a user-specified astropy 2D model.
The model parameters for each source are taken from the input
params_table
table.- Parameters:
- shape2-tuple of int
The shape of the output image.
- model2D
astropy.modeling.Model
The 2D model to be used to render the sources. The model must be two-dimensional where it accepts 2 inputs (i.e., (x, y)) and has 1 output. The model must have parameters for the x and y positions of the sources. Typically, these parameters are named ‘x_0’ and ‘y_0’, but the parameter names can be specified using the
x_name
andy_name
keywords.- params_table
Table
A table containing the model parameters for each source. Each row of the table corresponds to a source whose model parameters are defined by the column names, which must match the model parameter names. The table must contain columns for the x and y positions of the sources. The column names for the x and y positions can be specified using the
x_name
andy_name
keywords. Model parameters not defined in the table will be set to themodel
default value. To attach units to model parameters,params_table
must be input as aQTable
.If the table contains a column named ‘model_shape’, then the values in that column will be used to override the
model_shape
keyword and modelbounding_box
for each source. This can be used to render each source with a different shape.If the table contains a column named ‘local_bkg’, then the per-pixel local background values in that column will be used to added to each model source over the region defined by its
model_shape
. The ‘local_bkg’ column must have the same flux units as the output image (e.g., if the inputmodel
has ‘amplitude’ or ‘flux’ parameters with units). Including ‘local_bkg’ should be used with care, especially in crowded fields where themodel_shape
of sources overlap (see Notes below).Except for
model_shape
andlocal_bkg
column names, column names that do not match model parameters will be ignored.- model_shape2-tuple of int, int, or
None
, optional The shape around the (x, y) center of each source that will used to evaluate the
model
. Ifmodel_shape
is a scalar integer, then a square shape of sizemodel_shape
will be used. IfNone
, then the bounding box of the model will be used (which can optionally be scaled using thebbox_factor
keyword if the model supports it). This keyword must be specified if the model does not have abounding_box
attribute. If specified, this keyword overrides the modelbounding_box
attribute. To use a different shape for each source, include a column named'model_shape'
in theparams_table
. For that case, this keyword is ignored.- bbox_factor
None
or float, optional The multiplicative factor to pass to the model
bounding_box
method to determine the model shape. If the modelbounding_box
method does not accept afactor
keyword, then this keyword is ignored. IfNone
, the default model bounding box will be used. This keyword is ignored ifmodel_shape
is specified or if theparams_table
contains a'model_shape'
column. Note that some Photutils PSF models have abbox_factor
keyword that is be used to define the model bounding box. In that case, this keyword is ignored.- x_namestr, optional
The name of the
model
parameter that corresponds to the x position of the sources. This parameter must also be a column name inparams_table
.- y_namestr, optional
The name of the
model
parameter that corresponds to the y position of the sources. This parameter must also be a column name inparams_table
.- discretize_method{‘center’, ‘interp’, ‘oversample’, ‘integrate’}, optional
One of the following methods for discretizing the model on the pixel grid:
'center'
(default)Discretize model by taking the value at the center of the pixel bins. This method should be used for ePSF/PRF single or gridded models.
'interp'
Discretize model by bilinearly interpolating between the values at the corners of the pixel bins.
'oversample'
Discretize model by taking the average of model values in the pixel bins on an oversampled grid. Use the
discretize_oversample
keyword to set the integer oversampling factor.
'integrate'
Discretize model by integrating the model over the pixel bins using
scipy.integrate.quad
. This mode conserves the model integral on a subpixel scale, but it is extremely slow.
- discretize_oversampleint, optional
The integer oversampling factor used when
descretize_method='oversample'
. This keyword is ignored otherwise.- progress_barbool, optional
Whether to display a progress bar while adding the sources to the image. The progress bar requires that the tqdm optional dependency be installed. Note that the progress bar does not currently work in the Jupyter console due to limitations in
tqdm
.
- Returns:
- array2D
ndarray
The rendered image containing the model sources.
- array2D
Notes
The local background value around each source is optionally included using the
local_bkg
column in the inputparams_table
. This local background added to each source over itsmodel_shape
region. In regions where themodel_shape
of source overlap, the local background will be added multiple times. This is not an issue if the sources are well-separated, but for crowded fields, this option should be used with care.Examples
import matplotlib.pyplot as plt from astropy.modeling.models import Moffat2D from photutils.datasets import (make_model_image, make_random_models_table) model = Moffat2D() n_sources = 25 shape = (100, 100) param_ranges = {'amplitude': [100, 200], 'x_0': [0, shape[1]], 'y_0': [0, shape[0]], 'gamma': [1, 2], 'alpha': [1, 2]} params = make_random_models_table(n_sources, param_ranges, seed=0) model_shape = (15, 15) data = make_model_image(shape, model, params, model_shape=model_shape) plt.imshow(data, origin='lower') plt.tight_layout()
(
Source code
,png
,hires.png
,pdf
,svg
)import matplotlib.pyplot as plt import numpy as np from astropy.modeling.models import Gaussian2D from photutils.datasets import make_model_image, make_model_params model = Gaussian2D() shape = (500, 500) n_sources = 100 params = make_model_params(shape, n_sources, x_name='x_mean', y_name='y_mean', min_separation=25, amplitude=(100, 500), x_stddev=(1, 3), y_stddev=(1, 3), theta=(0, np.pi)) model_shape = (25, 25) data = make_model_image(shape, model, params, model_shape=model_shape, x_name='x_mean', y_name='y_mean') plt.imshow(data, origin='lower') plt.tight_layout()
(
Source code
,png
,hires.png
,pdf
,svg
)