18 from __future__
import division, print_function
25 import _SourceXtractorPy
as cpp
30 warnings.warn(
'Could not import pyston: running outside sourcextractor?', ImportWarning)
31 from astropy
import units
as u
32 from astropy.coordinates
import SkyCoord
34 from .measurement_images
import MeasurementGroup
44 Limit, and normalize, the range of values for a model fitting parameter.
49 limits : a tuple (min, max), or a callable that receives a source, and returns a tuple (min, max)
55 Normalized to engine space using a sigmoid function
59 engine = \ln \frac{world - min}{max-world} \\
60 world = min + \frac{max - min}{1 + e^{engine}}
63 Normalized to engine space using an exponential sigmoid function
67 engine = \ln \left( \frac{\ln(world/min)}{\ln(max /world)} \right) \\
68 world = min * e^\frac{ \ln(max / min) }{ (1 + e^{-engine}) }
77 self.
__callable = limits
if hasattr(limits,
'__call__')
else lambda v, o: limits
85 o : object being fitted
89 The minimum acceptable value for the range
98 o : object being fitted
102 The maximum acceptable value for the range
119 Human readable representation for the object
122 if hasattr(self.
__limits,
'__call__'):
127 RangeType.LINEAR:
'LIN',
128 RangeType.EXPONENTIAL:
'EXP'
130 res +=
',{}]'.format(type_str[self.
__type])
136 Unbounded, but normalize, value of a model fitting parameter
140 normalization_factor: float, or a callable that receives the initial value parameter value and a source,
142 The world value which will be normalized to 1 in engine coordinates
150 if hasattr(normalization_factor,
'__call__'):
153 self.
__callable =
lambda v, o: normalization_factor
160 o : object being fitted
165 The world value which will be normalized to 1 in engine coordinates
174 Human readable representation for the object
187 Base class for all model fitting parameter types.
188 Can not be used directly.
196 Human readable representation for the object
198 return '(ID:{})'.format(self.id)
203 A parameter with a single value that remains constant. It will not be fitted.
207 value : float, or callable that receives a source and returns a float
208 Value for this parameter
215 ParameterBase.__init__(self)
217 self.
__callable = value
if hasattr(value,
'__call__')
else lambda o: value
223 o : object being fitted
228 Value of the constant parameter
237 Human readable representation for the object
239 res = ParameterBase.__str__(self)[:-1] +
', value:'
240 if hasattr(self.
__value,
'__call__'):
249 A parameter that will be fitted by the model fitting engine.
253 init_value : float or callable that receives a source, and returns a float
254 Initial value for the parameter.
255 range : instance of Range or Unbounded
256 Defines if this parameter is unbounded or bounded, and how.
265 >>> sersic = FreeParameter(2.0, Range((1.0, 7.0), RangeType.LINEAR))
272 ParameterBase.__init__(self)
274 self.
__init_call = init_value
if hasattr(init_value,
'__call__')
else lambda o: init_value
281 o : object being fitted
286 Initial value for the free parameter
303 Human readable representation for the object
305 res = ParameterBase.__str__(self)[:-1] +
', init:'
311 res +=
', range:' + str(self.
__range)
317 A DependentParameter is not fitted by itself, but its value is derived from another Parameters, whatever their type:
318 FreeParameter, ConstantParameter, or other DependentParameter
323 A callable that will be called with all the parameters specified in this constructor each time a new
324 evaluation is needed.
325 params : list of ParameterBase
326 List of parameters on which this DependentParameter depends.
330 >>> flux = get_flux_parameter()
331 >>> mag = DependentParameter(lambda f: -2.5 * np.log10(f) + args.mag_zeropoint, flux)
332 >>> add_output_column('mf_mag_' + band, mag)
339 ParameterBase.__init__(self)
346 Convenience function for the position parameter X and Y.
351 X coordinate, starting at the X coordinate of the centroid and linearly limited to X +/- the object radius.
353 Y coordinate, starting at the Y coordinate of the centroid and linearly limited to Y +/- the object radius.
356 X and Y are fitted on the detection image X and Y coordinates. Internally, these are translated to measurement
357 images using the WCS headers.
359 param_range =
Range(
lambda v, o: (v - o.radius, v + o.radius), RangeType.LINEAR)
368 Possible flux types to use as initial value for the flux parameter.
369 Right now, only isophotal is supported.
376 Convenience function for the flux parameter.
381 One of the values defined in FluxParameterType
383 Scaling of the initial flux. Defaults to 1.
388 Flux parameter, starting at the flux defined by `type`, and limited to +/- 1e3 times the initial value.
391 FluxParameterType.ISO:
'isophotal_flux'
393 return FreeParameter(
lambda o: getattr(o, attr_map[type]) * scale,
394 Range(
lambda v, o: (v * 1E-3, v * 1E3), RangeType.EXPONENTIAL))
399 Model a Gaussian prior on a given parameter.
403 param : ParameterBase
404 Model fitting parameter
405 value : float or callable that receives a source and returns a float
407 sigma : float or callable that receives a source and returns a float
408 Standard deviation of the Gaussian
415 cpp.Id.__init__(self)
417 self.
value = value
if hasattr(value,
'__call__')
else lambda o: value
418 self.
sigma = sigma
if hasattr(sigma,
'__call__')
else lambda o: sigma
423 Base class for all models.
432 Base class for positioned models with a flux. It can not be used directly.
436 x_coord : ParameterBase or float
437 X coordinate (in the detection image)
438 y_coord : ParameterBase or float
439 Y coordinate (in the detection image)
440 flux : ParameterBase or float
448 ModelBase.__init__(self)
460 Models a source as a point, spread by the PSF.
464 x_coord : ParameterBase or float
465 X coordinate (in the detection image)
466 y_coord : ParameterBase or float
467 Y coordinate (in the detection image)
468 flux : ParameterBase or float
476 CoordinateModelBase.__init__(self, x_coord, y_coord, flux)
480 Return a human readable representation of the model.
485 If True, include information about the parameters.
492 return 'PointSource[x_coord={}, y_coord={}, flux={}]'.format(
495 return 'PointSource[x_coord={}, y_coord={}, flux={}]'.format(
496 self.x_coord.id, self.y_coord.id, self.flux.id)
501 A model that is constant through all the image.
505 value: ParameterBase or float
506 Value to add to the value of all pixels from the model.
513 ModelBase.__init__(self)
518 Return a human readable representation of the model.
523 If True, include information about the parameters.
530 return 'ConstantModel[value={}]'.format(self.
value)
532 return 'ConstantModel[value={}]'.format(self.value.id)
540 Base class for the Sersic, Exponential and de Vaucouleurs models. It can not be used directly.
544 x_coord : ParameterBase or float
545 X coordinate (in the detection image)
546 y_coord : ParameterBase or float
547 Y coordinate (in the detection image)
548 flux : ParameterBase or float
550 effective_radius : ParameterBase or float
551 Ellipse semi-major axis, in pixels on the detection image.
552 aspect_ratio : ParameterBase or float
554 angle : ParameterBase or float
555 Ellipse rotation, in radians
558 def __init__(self, x_coord, y_coord, flux, effective_radius, aspect_ratio, angle):
562 CoordinateModelBase.__init__(self, x_coord, y_coord, flux)
576 Model a source with a Sersic profile.
580 x_coord : ParameterBase or float
581 X coordinate (in the detection image)
582 y_coord : ParameterBase or float
583 Y coordinate (in the detection image)
584 flux : ParameterBase or float
586 effective_radius : ParameterBase or float
587 Ellipse semi-major axis, in pixels on the detection image.
588 aspect_ratio : ParameterBase or float
590 angle : ParameterBase or float
591 Ellipse rotation, in radians
592 n : ParameterBase or float
596 def __init__(self, x_coord, y_coord, flux, effective_radius, aspect_ratio, angle, n):
600 SersicModelBase.__init__(self, x_coord, y_coord, flux, effective_radius, aspect_ratio,
606 Return a human readable representation of the model.
611 If True, include information about the parameters.
618 return 'Sersic[x_coord={}, y_coord={}, flux={}, effective_radius={}, aspect_ratio={}, angle={}, n={}]'.format(
622 return 'Sersic[x_coord={}, y_coord={}, flux={}, effective_radius={}, aspect_ratio={}, angle={}, n={}]'.format(
623 self.x_coord.id, self.y_coord.id, self.flux.id, self.effective_radius.id,
624 self.aspect_ratio.id, self.angle.id, self.n.id)
627 return SersicModelBase.get_params(self) + [self.
n]
632 Model a source with an exponential profile (Sersic model with an index of 1)
636 x_coord : ParameterBase or float
637 X coordinate (in the detection image)
638 y_coord : ParameterBase or float
639 Y coordinate (in the detection image)
640 flux : ParameterBase or float
642 effective_radius : ParameterBase or float
643 Ellipse semi-major axis, in pixels on the detection image.
644 aspect_ratio : ParameterBase or float
646 angle : ParameterBase or float
647 Ellipse rotation, in radians
650 def __init__(self, x_coord, y_coord, flux, effective_radius, aspect_ratio, angle):
654 SersicModelBase.__init__(self, x_coord, y_coord, flux, effective_radius, aspect_ratio,
659 Return a human readable representation of the model.
664 If True, include information about the parameters.
671 return 'Exponential[x_coord={}, y_coord={}, flux={}, effective_radius={}, aspect_ratio={}, angle={}]'.format(
675 return 'Exponential[x_coord={}, y_coord={}, flux={}, effective_radius={}, aspect_ratio={}, angle={}]'.format(
676 self.x_coord.id, self.y_coord.id, self.flux.id, self.effective_radius.id,
677 self.aspect_ratio.id, self.angle.id)
682 Model a source with a De Vaucouleurs profile (Sersic model with an index of 4)
686 x_coord : ParameterBase or float
687 X coordinate (in the detection image)
688 y_coord : ParameterBase or float
689 Y coordinate (in the detection image)
690 flux : ParameterBase or float
692 effective_radius : ParameterBase or float
693 Ellipse semi-major axis, in pixels on the detection image.
694 aspect_ratio : ParameterBase or float
696 angle : ParameterBase or float
697 Ellipse rotation, in radians
700 def __init__(self, x_coord, y_coord, flux, effective_radius, aspect_ratio, angle):
704 SersicModelBase.__init__(self, x_coord, y_coord, flux, effective_radius, aspect_ratio,
709 Return a human readable representation of the model.
714 If True, include information about the parameters.
721 return 'DeVaucouleurs[x_coord={}, y_coord={}, flux={}, effective_radius={}, aspect_ratio={}, angle={}]'.format(
725 return 'DeVaucouleurs[x_coord={}, y_coord={}, flux={}, effective_radius={}, aspect_ratio={}, angle={}]'.format(
726 self.x_coord.id, self.y_coord.id, self.flux.id, self.effective_radius.id,
727 self.aspect_ratio.id, self.angle.id)
732 ComputeGraphModel model
736 models: string or list of strings, corresponding to path to Onnx format models,
737 (specifying more than one allows using the most efficient model based on render size.)
738 x_coord : ParameterBase or float
739 X coordinate (in the detection image)
740 y_coord : ParameterBase or float
741 Y coordinate (in the detection image)
742 flux : ParameterBase or float
744 params : Dictionary of String and ParameterBase or float
745 Dictionary of custom parameters for the ONNX model
748 def __init__(self, models, x_coord, y_coord, flux, params={}):
752 CoordinateModelBase.__init__(self, x_coord, y_coord, flux)
754 ratio_name =
"_aspect_ratio"
755 angle_name =
"_angle"
756 scale_name =
"_scale"
758 for k
in params.keys():
759 if not isinstance(params[k], ParameterBase):
762 aspect_ratio = params[ratio_name]
if ratio_name
in params.keys()
else 1.0
763 angle = params[angle_name]
if angle_name
in params.keys()
else 0.0
764 scale = params[scale_name]
if scale_name
in params.keys()
else 1.0
772 params.pop(ratio_name,
None)
773 params.pop(angle_name,
None)
774 params.pop(scale_name,
None)
777 self.
models = models
if isinstance(models, list)
else [models]
781 Return a human readable representation of the model.
786 If True, include information about the parameters.
793 return 'ComputeGraph[x_coord={}, y_coord={}, flux={}, effective_radius={}, aspect_ratio={}, angle={}]'.format(
797 return 'ComputeGraph[x_coord={}, y_coord={}, flux={}, effective_radius={}, aspect_ratio={}, angle={}]'.format(
798 self.x_coord.id, self.y_coord.id, self.flux.id, self.effective_radius.id,
799 self.aspect_ratio.id, self.angle.id)
803 list(self.params.values()))
808 Coordinates in right ascension and declination
828 Transform an (X, Y) in pixel coordinates on the detection image to (RA, DEC) in world coordinates.
839 wc_alpha = pyston.image_to_world_alpha(x - 1, y - 1)
840 wc_delta = pyston.image_to_world_delta(x - 1, y - 1)
846 Transform an (X, Y) in pixel coordinates on the detection image to astropy SkyCoord.
858 sky_coord = SkyCoord(ra=coord.ra * u.degree, dec=coord.dec * u.degree)
864 Transform a radius in pixels on the detection image to a radius in sky coordinates.
881 Get the separation angle in sky coordinates for two points defined in pixels on the detection image.
892 Separation in degrees
896 return c1.separation(c2).degree
901 Get the position angle in sky coordinates for two points defined in pixels on the detection image.
912 Position angle in degrees, normalized to -/+ 90
916 angle = c1.position_angle(c2).degree
919 return (angle + 90.0) % 180.0 - 90.0
924 Convenience function for generating two dependent parameter with world (alpha, delta) coordinates
925 from image (X, Y) coordinates.
934 ra : DependentParameter
935 dec : DependentParameter
943 >>> x, y = get_pos_parameters()
944 >>> ra, dec = get_world_position_parameters(x, y)
945 >>> add_output_column('mf_ra', ra)
946 >>> add_output_column('mf_dec', dec)
955 Convenience function for generating five dependent parameters, in world coordinates, for the position
956 and shape of a model.
962 radius : ParameterBase
963 angle : ParameterBase
964 ratio : ParameterBase
968 ra : DependentParameter
970 dec : DependentParameter
972 rad : DependentParameter
974 angle : DependentParameter
976 ratio : DependentParameter
977 Aspect ratio. It has to be recomputed as the axis of the ellipse may have different ratios
978 in image coordinates than in world coordinates
982 >>> flux = get_flux_parameter()
983 >>> x, y = get_pos_parameters()
984 >>> radius = FreeParameter(lambda o: o.radius, Range(lambda v, o: (.01 * v, 100 * v), RangeType.EXPONENTIAL))
985 >>> angle = FreeParameter(lambda o: o.angle, Range((-np.pi, np.pi), RangeType.LINEAR))
986 >>> ratio = FreeParameter(1, Range((0, 10), RangeType.LINEAR))
987 >>> add_model(group, ExponentialModel(x, y, flux, radius, ratio, angle))
988 >>> ra, dec, wc_rad, wc_angle, wc_ratio = get_world_parameters(x, y, radius, angle, ratio)
989 >>> add_output_column('mf_world_angle', wc_angle)
994 def get_major_axis(x, y, radius, angle, ratio):
996 x2 = x + math.cos(angle) * radius
997 y2 = y + math.sin(angle) * radius
999 x2 = x + math.sin(angle) * radius * ratio
1000 y2 = y + math.cos(angle) * radius * ratio
1004 def get_minor_axis(x, y, radius, angle, ratio):
1006 x2 = x + math.sin(angle) * radius * ratio
1007 y2 = y + math.cos(angle) * radius * ratio
1009 x2 = x + math.cos(angle) * radius
1010 y2 = y + math.sin(angle) * radius
1014 def wc_rad_func(x, y, radius, angle, ratio):
1015 x2, y2 = get_major_axis(x, y, radius, angle, ratio)
1018 def wc_angle_func(x, y, radius, angle, ratio):
1019 x2, y2 = get_major_axis(x, y, radius, angle, ratio)
1022 def wc_ratio_func(x, y, radius, angle, ratio):
1023 minor_x, minor_y = get_minor_axis(x, y, radius, angle, ratio)
1026 major_x, major_y = get_major_axis(x, y, radius, angle, ratio)
1029 return minor_angle / major_angle
1035 return ra, dec, wc_rad, wc_angle, wc_ratio
1051 self.
params_dict = {
"max_iterations": 200,
"modified_chi_squared_scale": 10,
"engine":
"",
1052 "use_iterative_fitting":
True,
"meta_iterations": 5,
1053 "deblend_factor": 0.95,
"meta_iteration_stop": 0.0001}
1057 if isinstance(x, tuple):
1070 if isinstance(attr, ConstantParameter):
1072 elif isinstance(attr, FreeParameter):
1074 elif isinstance(attr, DependentParameter):
1076 for param
in attr.params:
1080 for param
in model.get_params():
1084 if isinstance(model, ConstantModel):
1086 elif isinstance(model, PointSourceModel):
1088 elif isinstance(model, SersicModel):
1090 elif isinstance(model, ExponentialModel):
1092 elif isinstance(model, DeVaucouleursModel):
1094 elif isinstance(model, ComputeGraphModel):
1097 raise TypeError(
'Unknown model type {}'.format(type(model)))
1101 Add a model to be fitted to the given group.
1105 group : MeasurementGroup
1108 if not isinstance(group, MeasurementGroup):
1110 'Models can only be added on MeasurementGroup, got {}'.format(type(group)))
1111 if not hasattr(group,
'models'):
1113 group.models.append(model)
1120 Add a prior to the given parameter.
1124 param : ParameterBase
1125 value : float or callable that receives a source and returns a float
1126 Mean of the Gaussian
1127 sigma : float or callable that receives a source and returns a float
1128 Standard deviation of the Gaussian
1130 prior =
Prior(param, value, sigma)
1137 Print a human-readable representation of the configured model fitting parameters.
1142 Where to print the representation. Defaults to sys.stderr
1145 print(
'Constant parameters:', file=file)
1146 for n, param
in self.constant_parameter_dict.items():
1147 print(
' {}: {}'.format(n, param), file=file)
1149 print(
'Free parameters:', file=file)
1150 for n, param
in self.free_parameter_dict.items():
1151 print(
' {}: {}'.format(n, param), file=file)
1153 print(
'Dependent parameters:', file=file)
1154 for n, param
in self.dependent_parameter_dict.items():
1155 print(
' {}: {}'.format(n, param), file=file)
1162 Max number of iterations for the model fitting.
1171 Sets u0, as used by the modified chi squared residual comparator, a function that reduces the effect of large
1173 Refer to the SourceXtractor++ documentation for a better explanation of how residuals are computed and how
1174 this value affects the model fitting.
1176 self.
params_dict[
"modified_chi_squared_scale"] = scale
1183 Minimization engine for the model fitting : levmar or gsl
1191 use_iterative_fitting : boolean
1192 use iterative model fitting or legacy
1194 self.
params_dict[
"use_iterative_fitting"] = use_iterative_fitting
1200 meta_iterations : int
1201 number of meta iterations on the whole group (when using iterative model fitting)
1203 self.
params_dict[
"meta_iterations"] = meta_iterations
1211 self.
params_dict[
"deblend_factor"] = deblend_factor
1219 self.
params_dict[
"meta_iteration_stop"] = meta_iteration_stop
1224 Print a human-readable representation of the configured models.
1228 group : MeasurementGroup
1229 Print the models for this group.
1231 If True, print also the parameters that belong to the model
1233 Prefix each line with this string. Used internally for indentation.
1235 Where to print the representation. Defaults to sys.stderr
1237 if hasattr(group,
'models')
and group.models:
1238 print(
'{}Models:'.format(prefix), file=file)
1239 for m
in group.models:
1240 print(
'{} {}'.format(prefix, m.to_string(show_params)), file=file)
1242 if isinstance(x, tuple):
1243 print(
'{}{}:'.format(prefix, x[0]), file=file)