astropy.modeling provides a framework for representing models and performing model evaluation and fitting. It currently supports 1-D and 2-D models and fitting with parameter constraints.
It is designed to be easily extensible and flexible. Models do not reference fitting algorithms explicitly and new fitting algorithms may be added without changing the existing models (though not all models can be used with all fitting algorithms due to constraints such as model linearity).
The goal is to eventually provide a rich toolset of models and fitters such that most users will not need to define new model classes, nor special purpose fitting routines (while making it reasonably easy to do when necessary).
Warning
astropy.modeling is currently a work-in-progress, and thus it is likely there will be significant API changes in later versions of Astropy. If you have specific ideas for how it might be improved, feel free to let us know on the astropy-dev mailing list or at http://feedback.astropy.org
The examples here use the predefined models and assume the following modules have been imported:
>>> import numpy as np
>>> from astropy.modeling import models, fitting
The astropy.modeling package defines a number of models that are collected under a single namespace as astropy.modeling.models. Models behave like parametrized functions:
>>> from astropy.modeling import models
>>> g = models.Gaussian1D(amplitude=1.2, mean=0.9, stddev=0.5)
>>> print(g)
Model: Gaussian1D
Inputs: 1
Outputs: 1
Model set size: 1
Parameters:
amplitude mean stddev
--------- ---- ------
1.2 0.9 0.5
Model parameters can be accessed as attributes:
>>> g.amplitude
Parameter('amplitude', value=1.2)
>>> g.mean
Parameter('mean', value=0.9)
>>> g.stddev
Parameter('stddev', value=0.5)
and can also be updated via those attributes:
>>> g.amplitude = 0.8
>>> g.amplitude
Parameter('amplitude', value=0.8)
Models can be evaluated by calling them as functions:
>>> g(0.1)
0.22242984036255528
>>> g(np.linspace(0.5, 1.5, 7))
array([ 0.58091923, 0.71746405, 0.7929204 , 0.78415894, 0.69394278,
0.54952605, 0.3894018 ])
As the above example demonstrates, in general most models evaluate array-like inputs according to the standard Numpy broadcasting rules for arrays.
Models can therefore already be useful to evaluate common functions, independently of the fitting features of the package.
In this section, we look at a simple example of fitting a Gaussian to a simulated dataset. We use the Gaussian1D and Trapezoid1D models and the LevMarLSQFitter fitter to fit the data:
import numpy as np
from astropy.modeling import models, fitting
# Generate fake data
np.random.seed(0)
x = np.linspace(-5., 5., 200)
y = 3 * np.exp(-0.5 * (x - 1.3)**2 / 0.8**2)
y += np.random.normal(0., 0.2, x.shape)
# Fit the data using a box model
t_init = models.Trapezoid1D(amplitude=1., x_0=0., width=1., slope=0.5)
fit_t = fitting.LevMarLSQFitter()
t = fit_t(t_init, x, y)
# Fit the data using a Gaussian
g_init = models.Gaussian1D(amplitude=1., mean=0, stddev=1.)
fit_g = fitting.LevMarLSQFitter()
g = fit_g(g_init, x, y)
# Plot the data with the best-fit model
plt.figure(figsize=(8,5))
plt.plot(x, y, 'ko')
plt.plot(x, t(x), 'b-', lw=2, label='Trapezoid')
plt.plot(x, g(x), 'r-', lw=2, label='Gaussian')
plt.xlabel('Position')
plt.ylabel('Flux')
plt.legend(loc=2)
As shown above, once instantiated, the fitter class can be used as a function that takes the initial model (t_init or g_init) and the data values (x and y), and returns a fitted model (t or g).
Similarly to the 1-D example, we can create a simulated 2-D data dataset, and fit a polynomial model to it. This could be used for example to fit the background in an image.
import numpy as np
from astropy.modeling import models, fitting
# Generate fake data
np.random.seed(0)
y, x = np.mgrid[:128, :128]
z = 2. * x ** 2 - 0.5 * x ** 2 + 1.5 * x * y - 1.
z += np.random.normal(0., 0.1, z.shape) * 50000.
# Fit the data using astropy.modeling
p_init = models.Polynomial2D(degree=2)
fit_p = fitting.LevMarLSQFitter()
p = fit_p(p_init, x, y, z)
# Plot the data with the best-fit model
plt.figure(figsize=(8,2.5))
plt.subplot(1,3,1)
plt.imshow(z, interpolation='nearest', vmin=-1e4, vmax=5e4)
plt.title("Data")
plt.subplot(1,3,2)
plt.imshow(p(x, y), interpolation='nearest', vmin=-1e4, vmax=5e4)
plt.title("Model")
plt.subplot(1,3,3)
plt.imshow(z - p(x, y), interpolation='nearest', vmin=-1e4, vmax=5e4)
plt.title("Residual")
A list of models is provided in the Reference/API section. The fitting framework includes many useful features that are not demonstrated here, such as weighting of datapoints, fixing or linking parameters, and placing lower or upper limits on parameters. For more information on these, take a look at the Fitting Models to Data documentation.
In some cases it is necessary to describe many models of the same type but with different parameter values. This could be done simply by instantiating as many instances of a Model as are needed. But that can be inefficient for a large number of models. To that end, all model classes in astropy.modeling can also be used to represent a model set which is a collection of models of the same type, but with different values for their parameters.
To instantiate a model set, use argument n_models=N where N is the number of models in the set when constructing the model. The value of each parameter must be a list or array of length N, such that each item in the array corresponds to one model in the set:
>>> g = models.Gaussian1D(amplitude=[1, 2], mean=[0, 0],
... stddev=[0.1, 0.2], n_models=2)
>>> print(g)
Model: Gaussian1D
Inputs: 1
Outputs: 1
Model set size: 2
Parameters:
amplitude mean stddev
--------- ---- ------
1.0 0.0 0.1
2.0 0.0 0.2
This is equivalent to two Gaussians with the parameters amplitude=1, mean=0, stddev=0.1 and amplitude=2, mean=0, stddev=0.2 respectively. When printing the model the parameter values are displayed as a table, with each row corresponding to a single model in the set.
The number of models in a model set can be determined using the len builtin:
>>> len(g)
2
Single models have a length of 1, and are not considered a model set as such.
When evaluating a model set, by default the input must be the same length as the number of models, with one input per model:
>>> g([0, 0.1])
array([ 1. , 1.76499381])
The result is an array with one result per model in the set. It is also possible to broadcast a single value to all models in the set:
>>> g(0)
array([ 1., 2.])
Model sets are used primarily for fitting, allowing a large number of models of the same type to be fitted simultaneously (and independently from each other) to some large set of inputs. For example, fitting a polynomial to the time response of each pixel in a data cube. This can greatly speed up the fitting process, especially for linear models.
This subpackage provides a framework for representing models and performing model evaluation and fitting. It supports 1D and 2D models and fitting with parameter constraints. It has some predefined models and fitting routines.
format_input(func) | Wraps a model’s __call__ method so that the input arrays are converted into the appropriate shape given the model’s parameter dimensions. |
Fittable1DModel(*args, **kwargs) | Base class for one-dimensional fittable models. |
Fittable2DModel(*args, **kwargs) | Base class for one-dimensional fittable models. |
FittableModel(*args, **kwargs) | Base class for models that can be fitted using the built-in fitting algorithms. |
InputParameterError | Used for incorrect input parameter values and definitions. |
LabeledInput(data, labels) | Used by SerialCompositeModel and SummedCompositeModel to choose input data using labels. |
Model(*args, **kwargs) | Base class for all models. |
ModelDefinitionError | Used for incorrect models definitions |
Parameter([name, description, default, ...]) | Wraps individual parameters. |
SerialCompositeModel(transforms[, inmap, ...]) | Composite model that evaluates models in series. |
SummedCompositeModel(transforms[, inmap, outmap]) | Composite model that evaluates models in parallel. |
digraph inheritance7fcdcb79b9 { rankdir=LR; size="8.0, 12.0"; "Fittable1DModel" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.Fittable1DModel.html#astropy.modeling.Fittable1DModel",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Base class for one-dimensional fittable models.",height=0.25,shape=box,fontsize=10]; "FittableModel" -> "Fittable1DModel" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Fittable2DModel" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.Fittable2DModel.html#astropy.modeling.Fittable2DModel",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Base class for one-dimensional fittable models.",height=0.25,shape=box,fontsize=10]; "FittableModel" -> "Fittable2DModel" [arrowsize=0.5,style="setlinewidth(0.5)"]; "FittableModel" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.FittableModel.html#astropy.modeling.FittableModel",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Base class for models that can be fitted using the built-in fitting",height=0.25,shape=box,fontsize=10]; "Model" -> "FittableModel" [arrowsize=0.5,style="setlinewidth(0.5)"]; "InputParameterError" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.InputParameterError.html#astropy.modeling.InputParameterError",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Used for incorrect input parameter values and definitions.",height=0.25,shape=box,fontsize=10]; "LabeledInput" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.LabeledInput.html#astropy.modeling.LabeledInput",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Used by `SerialCompositeModel` and `SummedCompositeModel` to choose input",height=0.25,shape=box,fontsize=10]; "Model" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.Model.html#astropy.modeling.Model",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Base class for all models.",height=0.25,shape=box,fontsize=10]; "ModelDefinitionError" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.ModelDefinitionError.html#astropy.modeling.ModelDefinitionError",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Used for incorrect models definitions",height=0.25,shape=box,fontsize=10]; "Parameter" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.Parameter.html#astropy.modeling.Parameter",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Wraps individual parameters.",height=0.25,shape=box,fontsize=10]; "SerialCompositeModel" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.SerialCompositeModel.html#astropy.modeling.SerialCompositeModel",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Composite model that evaluates models in series.",height=0.25,shape=box,fontsize=10]; "_CompositeModel" -> "SerialCompositeModel" [arrowsize=0.5,style="setlinewidth(0.5)"]; "SummedCompositeModel" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.SummedCompositeModel.html#astropy.modeling.SummedCompositeModel",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Composite model that evaluates models in parallel.",height=0.25,shape=box,fontsize=10]; "_CompositeModel" -> "SummedCompositeModel" [arrowsize=0.5,style="setlinewidth(0.5)"]; "_CompositeModel" [shape=box,style="setlinewidth(0.5)",fontsize=10,fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",height=0.25]; "Model" -> "_CompositeModel" [arrowsize=0.5,style="setlinewidth(0.5)"]; }
This module implements classes (called Fitters) which combine optimization algorithms (typically from scipy.optimize) with statistic functions to perfom fitting. Fitters are implemented as callable classes. In addition to the data to fit, the __call__ method takes an instance of FittableModel as input, and returns a copy of the model with its parameters determined by the optimizer.
Optimization algorithms, called “optimizers” are implemented in optimizers and statistic functions are in statistic. The goal is to provide an easy to extend framework and allow users to easily create new fitters by combining statistics with optimizers.
There are two exceptions to the above scheme. LinearLSQFitter uses Numpy’s lstsq function. LevMarLSQFitter uses leastsq which combines optimization and statistic in one implementation.
LinearLSQFitter() | A class performing a linear least square fitting. |
LevMarLSQFitter() | Levenberg-Marquardt algorithm and least squares statistic. |
SLSQPLSQFitter() | SLSQP optimization algorithm and least squares statistic. |
SimplexLSQFitter() | Simplex algorithm and least squares statistic. |
JointFitter(models, jointparameters, initvals) | Fit models which share a parameter. |
Fitter(optimizer, statistic) | Base class for all fitters. |
digraph inheritance2d0d457931 { rankdir=LR; size="8.0, 12.0"; "Fitter" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.fitting.Fitter.html#astropy.modeling.fitting.Fitter",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Base class for all fitters.",height=0.25,shape=box,fontsize=10]; "JointFitter" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.fitting.JointFitter.html#astropy.modeling.fitting.JointFitter",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Fit models which share a parameter.",height=0.25,shape=box,fontsize=10]; "LevMarLSQFitter" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.fitting.LevMarLSQFitter.html#astropy.modeling.fitting.LevMarLSQFitter",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Levenberg-Marquardt algorithm and least squares statistic.",height=0.25,shape=box,fontsize=10]; "LinearLSQFitter" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.fitting.LinearLSQFitter.html#astropy.modeling.fitting.LinearLSQFitter",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="A class performing a linear least square fitting.",height=0.25,shape=box,fontsize=10]; "SLSQPLSQFitter" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.fitting.SLSQPLSQFitter.html#astropy.modeling.fitting.SLSQPLSQFitter",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="SLSQP optimization algorithm and least squares statistic.",height=0.25,shape=box,fontsize=10]; "Fitter" -> "SLSQPLSQFitter" [arrowsize=0.5,style="setlinewidth(0.5)"]; "SimplexLSQFitter" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.fitting.SimplexLSQFitter.html#astropy.modeling.fitting.SimplexLSQFitter",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Simplex algorithm and least squares statistic.",height=0.25,shape=box,fontsize=10]; "Fitter" -> "SimplexLSQFitter" [arrowsize=0.5,style="setlinewidth(0.5)"]; }
Mathematical models.
custom_model_1d(func[, func_fit_deriv]) | Create a one dimensional model from a user defined function. |
AiryDisk2D(amplitude, x_0, y_0, radius, **kwargs) | Two dimensional Airy disk model. |
Beta1D(amplitude, x_0, gamma, alpha, **kwargs) | One dimensional Beta model. |
Beta2D(amplitude, x_0, y_0, gamma, alpha, ...) | Two dimensional Beta model. |
Box1D(amplitude, x_0, width, **kwargs) | One dimensional Box model. |
Box2D(amplitude, x_0, y_0, x_width, y_width, ...) | Two dimensional Box model. |
Const1D(amplitude, **kwargs) | One dimensional Constant model. |
Const2D(amplitude, **kwargs) | Two dimensional Constant model. |
Disk2D(amplitude, x_0, y_0, R_0, **kwargs) | Two dimensional radial symmetric Disk model. |
Gaussian1D(amplitude, mean, stddev, **kwargs) | One dimensional Gaussian model. |
Gaussian2D(amplitude, x_mean, y_mean[, ...]) | Two dimensional Gaussian model. |
GaussianAbsorption1D(amplitude, mean, ...) | One dimensional Gaussian absorption line model. |
Linear1D(slope, intercept, **kwargs) | One dimensional Line model. |
Lorentz1D(amplitude, x_0, fwhm, **kwargs) | One dimensional Lorentzian model. |
MexicanHat1D(amplitude, x_0, sigma, **kwargs) | One dimensional Mexican Hat model. |
MexicanHat2D(amplitude, x_0, y_0, sigma, ...) | Two dimensional symmetric Mexican Hat model. |
Redshift(z, **kwargs) | One dimensional redshift model. |
Ring2D(amplitude, x_0, y_0, r_in[, width, r_out]) | Two dimensional radial symmetric Ring model. |
Scale(factors, **kwargs) | Multiply a model by a factor. |
Shift(offsets, **kwargs) | Shift a coordinate. |
Sine1D(amplitude, frequency, **kwargs) | One dimensional Sine model. |
Trapezoid1D(amplitude, x_0, width, slope, ...) | One dimensional Trapezoid model. |
TrapezoidDisk2D(amplitude, x_0, y_0, R_0, ...) | Two dimensional circular Trapezoid model. |
digraph inheritance5f754a30cf { rankdir=LR; size="8.0, 12.0"; "AiryDisk2D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.functional_models.AiryDisk2D.html#astropy.modeling.functional_models.AiryDisk2D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Two dimensional Airy disk model.",height=0.25,shape=box,fontsize=10]; "Fittable2DModel" -> "AiryDisk2D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Beta1D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.functional_models.Beta1D.html#astropy.modeling.functional_models.Beta1D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="One dimensional Beta model.",height=0.25,shape=box,fontsize=10]; "Fittable1DModel" -> "Beta1D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Beta2D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.functional_models.Beta2D.html#astropy.modeling.functional_models.Beta2D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Two dimensional Beta model.",height=0.25,shape=box,fontsize=10]; "Fittable2DModel" -> "Beta2D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Box1D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.functional_models.Box1D.html#astropy.modeling.functional_models.Box1D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="One dimensional Box model.",height=0.25,shape=box,fontsize=10]; "Fittable1DModel" -> "Box1D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Box2D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.functional_models.Box2D.html#astropy.modeling.functional_models.Box2D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Two dimensional Box model.",height=0.25,shape=box,fontsize=10]; "Fittable2DModel" -> "Box2D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Const1D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.functional_models.Const1D.html#astropy.modeling.functional_models.Const1D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="One dimensional Constant model.",height=0.25,shape=box,fontsize=10]; "Fittable1DModel" -> "Const1D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Const2D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.functional_models.Const2D.html#astropy.modeling.functional_models.Const2D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Two dimensional Constant model.",height=0.25,shape=box,fontsize=10]; "Fittable2DModel" -> "Const2D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Disk2D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.functional_models.Disk2D.html#astropy.modeling.functional_models.Disk2D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Two dimensional radial symmetric Disk model.",height=0.25,shape=box,fontsize=10]; "Fittable2DModel" -> "Disk2D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Fittable1DModel" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.Fittable1DModel.html#astropy.modeling.Fittable1DModel",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Base class for one-dimensional fittable models.",height=0.25,shape=box,fontsize=10]; "FittableModel" -> "Fittable1DModel" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Fittable2DModel" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.Fittable2DModel.html#astropy.modeling.Fittable2DModel",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Base class for one-dimensional fittable models.",height=0.25,shape=box,fontsize=10]; "FittableModel" -> "Fittable2DModel" [arrowsize=0.5,style="setlinewidth(0.5)"]; "FittableModel" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.FittableModel.html#astropy.modeling.FittableModel",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Base class for models that can be fitted using the built-in fitting",height=0.25,shape=box,fontsize=10]; "Model" -> "FittableModel" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Gaussian1D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.functional_models.Gaussian1D.html#astropy.modeling.functional_models.Gaussian1D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="One dimensional Gaussian model.",height=0.25,shape=box,fontsize=10]; "Fittable1DModel" -> "Gaussian1D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Gaussian2D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.functional_models.Gaussian2D.html#astropy.modeling.functional_models.Gaussian2D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Two dimensional Gaussian model.",height=0.25,shape=box,fontsize=10]; "Fittable2DModel" -> "Gaussian2D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "GaussianAbsorption1D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.functional_models.GaussianAbsorption1D.html#astropy.modeling.functional_models.GaussianAbsorption1D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="One dimensional Gaussian absorption line model.",height=0.25,shape=box,fontsize=10]; "Fittable1DModel" -> "GaussianAbsorption1D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Linear1D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.functional_models.Linear1D.html#astropy.modeling.functional_models.Linear1D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="One dimensional Line model.",height=0.25,shape=box,fontsize=10]; "Fittable1DModel" -> "Linear1D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Lorentz1D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.functional_models.Lorentz1D.html#astropy.modeling.functional_models.Lorentz1D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="One dimensional Lorentzian model.",height=0.25,shape=box,fontsize=10]; "Fittable1DModel" -> "Lorentz1D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "MexicanHat1D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.functional_models.MexicanHat1D.html#astropy.modeling.functional_models.MexicanHat1D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="One dimensional Mexican Hat model.",height=0.25,shape=box,fontsize=10]; "Fittable1DModel" -> "MexicanHat1D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "MexicanHat2D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.functional_models.MexicanHat2D.html#astropy.modeling.functional_models.MexicanHat2D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Two dimensional symmetric Mexican Hat model.",height=0.25,shape=box,fontsize=10]; "Fittable2DModel" -> "MexicanHat2D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Model" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.Model.html#astropy.modeling.Model",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Base class for all models.",height=0.25,shape=box,fontsize=10]; "Redshift" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.functional_models.Redshift.html#astropy.modeling.functional_models.Redshift",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="One dimensional redshift model.",height=0.25,shape=box,fontsize=10]; "Fittable1DModel" -> "Redshift" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Ring2D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.functional_models.Ring2D.html#astropy.modeling.functional_models.Ring2D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Two dimensional radial symmetric Ring model.",height=0.25,shape=box,fontsize=10]; "Fittable2DModel" -> "Ring2D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Scale" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.functional_models.Scale.html#astropy.modeling.functional_models.Scale",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Multiply a model by a factor.",height=0.25,shape=box,fontsize=10]; "Model" -> "Scale" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Shift" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.functional_models.Shift.html#astropy.modeling.functional_models.Shift",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Shift a coordinate.",height=0.25,shape=box,fontsize=10]; "Model" -> "Shift" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Sine1D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.functional_models.Sine1D.html#astropy.modeling.functional_models.Sine1D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="One dimensional Sine model.",height=0.25,shape=box,fontsize=10]; "Fittable1DModel" -> "Sine1D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Trapezoid1D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.functional_models.Trapezoid1D.html#astropy.modeling.functional_models.Trapezoid1D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="One dimensional Trapezoid model.",height=0.25,shape=box,fontsize=10]; "Fittable1DModel" -> "Trapezoid1D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "TrapezoidDisk2D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.functional_models.TrapezoidDisk2D.html#astropy.modeling.functional_models.TrapezoidDisk2D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Two dimensional circular Trapezoid model.",height=0.25,shape=box,fontsize=10]; "Fittable2DModel" -> "TrapezoidDisk2D" [arrowsize=0.5,style="setlinewidth(0.5)"]; }
Optimization algorithms used in fitting.
Optimization(opt_method) | Base class for optimizers. |
SLSQP() | Sequential Least Squares Programming optimization algorithm. |
Simplex() | Neald-Mead (downhill simplex) algorithm [1]. |
digraph inheritancee60c3de4ec { rankdir=LR; size="8.0, 12.0"; "Optimization" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.optimizers.Optimization.html#astropy.modeling.optimizers.Optimization",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Base class for optimizers.",height=0.25,shape=box,fontsize=10]; "SLSQP" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.optimizers.SLSQP.html#astropy.modeling.optimizers.SLSQP",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Sequential Least Squares Programming optimization algorithm.",height=0.25,shape=box,fontsize=10]; "Optimization" -> "SLSQP" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Simplex" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.optimizers.Simplex.html#astropy.modeling.optimizers.Simplex",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Neald-Mead (downhill simplex) algorithm [1].",height=0.25,shape=box,fontsize=10]; "Optimization" -> "Simplex" [arrowsize=0.5,style="setlinewidth(0.5)"]; }
Power law model variants
BrokenPowerLaw1D(amplitude, x_break, ...) | One dimensional power law model with a break. |
ExponentialCutoffPowerLaw1D(amplitude, x_0, ...) | One dimensional power law model with an exponential cutoff. |
LogParabola1D(amplitude, x_0, alpha, beta, ...) | One dimensional log parabola model (sometimes called curved power law). |
PowerLaw1D(amplitude, x_0, alpha, **constraints) | One dimensional power law model. |
digraph inheritanceb72993dc76 { rankdir=LR; size="8.0, 12.0"; "BrokenPowerLaw1D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.powerlaws.BrokenPowerLaw1D.html#astropy.modeling.powerlaws.BrokenPowerLaw1D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="One dimensional power law model with a break.",height=0.25,shape=box,fontsize=10]; "Fittable1DModel" -> "BrokenPowerLaw1D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "ExponentialCutoffPowerLaw1D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.powerlaws.ExponentialCutoffPowerLaw1D.html#astropy.modeling.powerlaws.ExponentialCutoffPowerLaw1D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="One dimensional power law model with an exponential cutoff.",height=0.25,shape=box,fontsize=10]; "Fittable1DModel" -> "ExponentialCutoffPowerLaw1D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Fittable1DModel" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.Fittable1DModel.html#astropy.modeling.Fittable1DModel",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Base class for one-dimensional fittable models.",height=0.25,shape=box,fontsize=10]; "FittableModel" -> "Fittable1DModel" [arrowsize=0.5,style="setlinewidth(0.5)"]; "FittableModel" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.FittableModel.html#astropy.modeling.FittableModel",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Base class for models that can be fitted using the built-in fitting",height=0.25,shape=box,fontsize=10]; "Model" -> "FittableModel" [arrowsize=0.5,style="setlinewidth(0.5)"]; "LogParabola1D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.powerlaws.LogParabola1D.html#astropy.modeling.powerlaws.LogParabola1D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="One dimensional log parabola model (sometimes called curved power law).",height=0.25,shape=box,fontsize=10]; "Fittable1DModel" -> "LogParabola1D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Model" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.Model.html#astropy.modeling.Model",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Base class for all models.",height=0.25,shape=box,fontsize=10]; "PowerLaw1D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.powerlaws.PowerLaw1D.html#astropy.modeling.powerlaws.PowerLaw1D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="One dimensional power law model.",height=0.25,shape=box,fontsize=10]; "Fittable1DModel" -> "PowerLaw1D" [arrowsize=0.5,style="setlinewidth(0.5)"]; }
This module contains predefined polynomial models.
Chebyshev1D(degree[, domain, window, ...]) | 1D Chebyshev polynomial of the 1st kind. |
Chebyshev2D(x_degree, y_degree[, x_domain, ...]) | 2D Chebyshev polynomial of the 1st kind. |
InverseSIP(ap_order, bp_order[, ap_coeff, ...]) | Inverse Simple Imaging Polynomial |
Legendre1D(degree[, domain, window, ...]) | 1D Legendre polynomial. |
Legendre2D(x_degree, y_degree[, x_domain, ...]) | Legendre 2D polynomial. |
Polynomial1D(degree[, domain, window, ...]) | 1D Polynomial model. |
Polynomial2D(degree[, x_domain, y_domain, ...]) | 2D Polynomial model. |
SIP(crpix, a_order, b_order[, a_coeff, ...]) | Simple Imaging Polynomial (SIP) model. |
OrthoPolynomialBase(x_degree, y_degree[, ...]) | This is a base class for the 2D Chebyshev and Legendre models. |
PolynomialModel(degree[, n_inputs, ...]) | Base class for polynomial models. |
digraph inheritancecc8b866a55 { rankdir=LR; size="8.0, 12.0"; "Chebyshev1D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.polynomial.Chebyshev1D.html#astropy.modeling.polynomial.Chebyshev1D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="1D Chebyshev polynomial of the 1st kind.",height=0.25,shape=box,fontsize=10]; "PolynomialModel" -> "Chebyshev1D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Chebyshev2D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.polynomial.Chebyshev2D.html#astropy.modeling.polynomial.Chebyshev2D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="2D Chebyshev polynomial of the 1st kind.",height=0.25,shape=box,fontsize=10]; "OrthoPolynomialBase" -> "Chebyshev2D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "FittableModel" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.FittableModel.html#astropy.modeling.FittableModel",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Base class for models that can be fitted using the built-in fitting",height=0.25,shape=box,fontsize=10]; "Model" -> "FittableModel" [arrowsize=0.5,style="setlinewidth(0.5)"]; "InverseSIP" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.polynomial.InverseSIP.html#astropy.modeling.polynomial.InverseSIP",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Inverse Simple Imaging Polynomial",height=0.25,shape=box,fontsize=10]; "Model" -> "InverseSIP" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Legendre1D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.polynomial.Legendre1D.html#astropy.modeling.polynomial.Legendre1D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="1D Legendre polynomial.",height=0.25,shape=box,fontsize=10]; "PolynomialModel" -> "Legendre1D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Legendre2D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.polynomial.Legendre2D.html#astropy.modeling.polynomial.Legendre2D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Legendre 2D polynomial.",height=0.25,shape=box,fontsize=10]; "OrthoPolynomialBase" -> "Legendre2D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Model" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.Model.html#astropy.modeling.Model",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Base class for all models.",height=0.25,shape=box,fontsize=10]; "OrthoPolynomialBase" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.polynomial.OrthoPolynomialBase.html#astropy.modeling.polynomial.OrthoPolynomialBase",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="This is a base class for the 2D Chebyshev and Legendre models.",height=0.25,shape=box,fontsize=10]; "PolynomialBase" -> "OrthoPolynomialBase" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Polynomial1D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.polynomial.Polynomial1D.html#astropy.modeling.polynomial.Polynomial1D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="1D Polynomial model.",height=0.25,shape=box,fontsize=10]; "PolynomialModel" -> "Polynomial1D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Polynomial2D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.polynomial.Polynomial2D.html#astropy.modeling.polynomial.Polynomial2D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="2D Polynomial model.",height=0.25,shape=box,fontsize=10]; "PolynomialModel" -> "Polynomial2D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "PolynomialBase" [style="setlinewidth(0.5)",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Base class for all polynomial-like models with an arbitrary number of",height=0.25,shape=box,fontsize=10]; "FittableModel" -> "PolynomialBase" [arrowsize=0.5,style="setlinewidth(0.5)"]; "PolynomialModel" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.polynomial.PolynomialModel.html#astropy.modeling.polynomial.PolynomialModel",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Base class for polynomial models.",height=0.25,shape=box,fontsize=10]; "PolynomialBase" -> "PolynomialModel" [arrowsize=0.5,style="setlinewidth(0.5)"]; "SIP" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.polynomial.SIP.html#astropy.modeling.polynomial.SIP",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Simple Imaging Polynomial (SIP) model.",height=0.25,shape=box,fontsize=10]; "Model" -> "SIP" [arrowsize=0.5,style="setlinewidth(0.5)"]; }
Implements projections–particularly sky projections defined in WCS Paper II [R19]
All angles are set and and displayed in degrees but internally computations are performed in radians.
Pix2Sky_AZP([mu, gamma]) | AZP : Zenital perspective projection - pixel to sky. |
Sky2Pix_AZP([mu, gamma]) | AZP : Zenital perspective projection - sky to pixel. |
Pix2Sky_CAR(*args, **kwargs) | CAR: Plate carree projection - pixel to sky. |
Sky2Pix_CAR(*args, **kwargs) | CAR: Plate carree projection - sky to pixel. |
Pix2Sky_CEA([lam]) | CEA : Cylindrical equal area projection - pixel to sky. |
Sky2Pix_CEA([lam]) | CEA: Cylindrical equal area projection - sky to pixel. |
Pix2Sky_CYP(mu, lam) | CYP : Cylindrical perspective - pixel to sky. |
Sky2Pix_CYP(mu, lam) | CYP : Cylindrical Perspective - sky to pixel. |
Pix2Sky_MER(*args, **kwargs) | MER: Mercator - pixel to sky. |
Sky2Pix_MER(*args, **kwargs) | MER: Mercator - sky to pixel. |
Pix2Sky_SIN(*args, **kwargs) | SIN : Slant orthographic projection - pixel to sky. |
Sky2Pix_SIN(*args, **kwargs) | SIN : Slant othographic projection - sky to pixel. |
Pix2Sky_STG(*args, **kwargs) | STG : Stereographic Projection - pixel to sky. |
Sky2Pix_STG(*args, **kwargs) | STG : Stereographic Projection - sky to pixel. |
Pix2Sky_TAN(*args, **kwargs) | TAN : Gnomonic projection - pixel to sky. |
Sky2Pix_TAN(*args, **kwargs) | TAN : Gnomonic Projection - sky to pixel. |
AffineTransformation2D([matrix, translation]) | Perform an affine transformation in 2 dimensions. |
digraph inheritance447b6afb60 { rankdir=LR; size="8.0, 12.0"; "AffineTransformation2D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.projections.AffineTransformation2D.html#astropy.modeling.projections.AffineTransformation2D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Perform an affine transformation in 2 dimensions.",height=0.25,shape=box,fontsize=10]; "Model" -> "AffineTransformation2D" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Cylindrical" [style="setlinewidth(0.5)",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Base class for Cylindrical projections.",height=0.25,shape=box,fontsize=10]; "Projection" -> "Cylindrical" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Model" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.Model.html#astropy.modeling.Model",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Base class for all models.",height=0.25,shape=box,fontsize=10]; "Pix2Sky_AZP" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.projections.Pix2Sky_AZP.html#astropy.modeling.projections.Pix2Sky_AZP",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="AZP : Zenital perspective projection - pixel to sky.",height=0.25,shape=box,fontsize=10]; "Zenithal" -> "Pix2Sky_AZP" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Pix2Sky_CAR" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.projections.Pix2Sky_CAR.html#astropy.modeling.projections.Pix2Sky_CAR",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="CAR: Plate carree projection - pixel to sky.",height=0.25,shape=box,fontsize=10]; "Cylindrical" -> "Pix2Sky_CAR" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Pix2Sky_CEA" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.projections.Pix2Sky_CEA.html#astropy.modeling.projections.Pix2Sky_CEA",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="CEA : Cylindrical equal area projection - pixel to sky.",height=0.25,shape=box,fontsize=10]; "Cylindrical" -> "Pix2Sky_CEA" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Pix2Sky_CYP" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.projections.Pix2Sky_CYP.html#astropy.modeling.projections.Pix2Sky_CYP",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="CYP : Cylindrical perspective - pixel to sky.",height=0.25,shape=box,fontsize=10]; "Cylindrical" -> "Pix2Sky_CYP" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Pix2Sky_MER" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.projections.Pix2Sky_MER.html#astropy.modeling.projections.Pix2Sky_MER",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="MER: Mercator - pixel to sky.",height=0.25,shape=box,fontsize=10]; "Cylindrical" -> "Pix2Sky_MER" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Pix2Sky_SIN" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.projections.Pix2Sky_SIN.html#astropy.modeling.projections.Pix2Sky_SIN",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="SIN : Slant orthographic projection - pixel to sky.",height=0.25,shape=box,fontsize=10]; "Zenithal" -> "Pix2Sky_SIN" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Pix2Sky_STG" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.projections.Pix2Sky_STG.html#astropy.modeling.projections.Pix2Sky_STG",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="STG : Stereographic Projection - pixel to sky.",height=0.25,shape=box,fontsize=10]; "Zenithal" -> "Pix2Sky_STG" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Pix2Sky_TAN" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.projections.Pix2Sky_TAN.html#astropy.modeling.projections.Pix2Sky_TAN",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="TAN : Gnomonic projection - pixel to sky.",height=0.25,shape=box,fontsize=10]; "Zenithal" -> "Pix2Sky_TAN" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Projection" [style="setlinewidth(0.5)",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Base class for all sky projections.",height=0.25,shape=box,fontsize=10]; "Model" -> "Projection" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Sky2Pix_AZP" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.projections.Sky2Pix_AZP.html#astropy.modeling.projections.Sky2Pix_AZP",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="AZP : Zenital perspective projection - sky to pixel.",height=0.25,shape=box,fontsize=10]; "Zenithal" -> "Sky2Pix_AZP" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Sky2Pix_CAR" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.projections.Sky2Pix_CAR.html#astropy.modeling.projections.Sky2Pix_CAR",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="CAR: Plate carree projection - sky to pixel.",height=0.25,shape=box,fontsize=10]; "Cylindrical" -> "Sky2Pix_CAR" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Sky2Pix_CEA" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.projections.Sky2Pix_CEA.html#astropy.modeling.projections.Sky2Pix_CEA",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="CEA: Cylindrical equal area projection - sky to pixel.",height=0.25,shape=box,fontsize=10]; "Cylindrical" -> "Sky2Pix_CEA" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Sky2Pix_CYP" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.projections.Sky2Pix_CYP.html#astropy.modeling.projections.Sky2Pix_CYP",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="CYP : Cylindrical Perspective - sky to pixel.",height=0.25,shape=box,fontsize=10]; "Cylindrical" -> "Sky2Pix_CYP" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Sky2Pix_MER" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.projections.Sky2Pix_MER.html#astropy.modeling.projections.Sky2Pix_MER",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="MER: Mercator - sky to pixel.",height=0.25,shape=box,fontsize=10]; "Cylindrical" -> "Sky2Pix_MER" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Sky2Pix_SIN" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.projections.Sky2Pix_SIN.html#astropy.modeling.projections.Sky2Pix_SIN",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="SIN : Slant othographic projection - sky to pixel.",height=0.25,shape=box,fontsize=10]; "Zenithal" -> "Sky2Pix_SIN" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Sky2Pix_STG" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.projections.Sky2Pix_STG.html#astropy.modeling.projections.Sky2Pix_STG",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="STG : Stereographic Projection - sky to pixel.",height=0.25,shape=box,fontsize=10]; "Zenithal" -> "Sky2Pix_STG" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Sky2Pix_TAN" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.projections.Sky2Pix_TAN.html#astropy.modeling.projections.Sky2Pix_TAN",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="TAN : Gnomonic Projection - sky to pixel.",height=0.25,shape=box,fontsize=10]; "Zenithal" -> "Sky2Pix_TAN" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Zenithal" [style="setlinewidth(0.5)",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Base class for all Zenithal projections.",height=0.25,shape=box,fontsize=10]; "Projection" -> "Zenithal" [arrowsize=0.5,style="setlinewidth(0.5)"]; }
Statistic functions used in fitting.
leastsquare(measured_vals, updated_model, ...) | Least square statistic with optional weights. |
Implements roations, including spherical rotations as defined in WCS Paper II [R20]
RotateNative2Celestial and RotateCelestial2Native follow the convention in WCS Paper II to rotate to/from a native sphere and the celestial sphere.
The user interface sets and displays angles in degrees but the values are stored internally in radians. This is managed through the parameter setters/getters.
RotateCelestial2Native(phi, theta, psi) | Transformation from Celestial to Native to Spherical Coordinates. |
RotateNative2Celestial(phi, theta, psi) | Transformation from Native to Celestial Spherical Coordinates. |
Rotation2D([angle]) | Perform a 2D rotation given an angle in degrees. |
digraph inheritance5a3416a143 { rankdir=LR; size="8.0, 12.0"; "EulerAngleRotation" [style="setlinewidth(0.5)",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Base class for Euler angle rotations.",height=0.25,shape=box,fontsize=10]; "Model" -> "EulerAngleRotation" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Model" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.Model.html#astropy.modeling.Model",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Base class for all models.",height=0.25,shape=box,fontsize=10]; "RotateCelestial2Native" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.rotations.RotateCelestial2Native.html#astropy.modeling.rotations.RotateCelestial2Native",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Transformation from Celestial to Native to Spherical Coordinates.",height=0.25,shape=box,fontsize=10]; "EulerAngleRotation" -> "RotateCelestial2Native" [arrowsize=0.5,style="setlinewidth(0.5)"]; "RotateNative2Celestial" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.rotations.RotateNative2Celestial.html#astropy.modeling.rotations.RotateNative2Celestial",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Transformation from Native to Celestial Spherical Coordinates.",height=0.25,shape=box,fontsize=10]; "EulerAngleRotation" -> "RotateNative2Celestial" [arrowsize=0.5,style="setlinewidth(0.5)"]; "Rotation2D" [style="setlinewidth(0.5)",URL="../api/astropy.modeling.rotations.Rotation2D.html#astropy.modeling.rotations.Rotation2D",fontname="Vera Sans, DejaVu Sans, Liberation Sans, Arial, Helvetica, sans",tooltip="Perform a 2D rotation given an angle in degrees.",height=0.25,shape=box,fontsize=10]; "Model" -> "Rotation2D" [arrowsize=0.5,style="setlinewidth(0.5)"]; }