SourceXtractorPlusPlus 1.0.3
SourceXtractor++, the next generation SExtractor
Loading...
Searching...
No Matches
model_fitting.py
Go to the documentation of this file.
1# -*- coding: utf-8 -*-
2
3# Copyright © 2019 Université de Genève, LMU Munich - Faculty of Physics, IAP-CNRS/Sorbonne Université
4#
5# This library is free software; you can redistribute it and/or modify it under
6# the terms of the GNU Lesser General Public License as published by the Free
7# Software Foundation; either version 3.0 of the License, or (at your option)
8# any later version.
9#
10# This library is distributed in the hope that it will be useful, but WITHOUT
11# ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
12# FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more
13# details.
14#
15# You should have received a copy of the GNU Lesser General Public License
16# along with this library; if not, write to the Free Software Foundation, Inc.,
17# 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
18from __future__ import division, print_function
19
20import math
21import sys
22import warnings
23from enum import Enum
24
25import _SourceXtractorPy as cpp
26
27try:
28 import pyston
29except ImportError:
30 warnings.warn('Could not import pyston: running outside sourcextractor?', ImportWarning)
31from astropy import units as u
32from astropy.coordinates import SkyCoord
33
34from .measurement_images import MeasurementGroup
35
36
37class RangeType(Enum):
38 LINEAR = 1
39 EXPONENTIAL = 2
40
41
42class Range(object):
43 r"""
44 Limit, and normalize, the range of values for a model fitting parameter.
45
46
47 Parameters
48 ----------
49 limits : a tuple (min, max), or a callable that receives a source, and returns a tuple (min, max)
50 type : RangeType
51
52 Notes
53 -----
54 RangeType.LINEAR
55 Normalized to engine space using a sigmoid function
56
57 .. math::
58
59 engine = \ln \frac{world - min}{max-world} \\
60 world = min + \frac{max - min}{1 + e^{engine}}
61
62 RangeType.EXPONENTIAL
63 Normalized to engine space using an exponential sigmoid function
64
65 .. math::
66
67 engine = \ln \left( \frac{\ln(world/min)}{\ln(max /world)} \right) \\
68 world = min * e^\frac{ \ln(max / min) }{ (1 + e^{-engine}) }
69
70 """
71
72 def __init__(self, limits, type):
73 """
74 Constructor.
75 """
76 self.__limits = limits
77 self.__callable = limits if hasattr(limits, '__call__') else lambda v, o: limits
78 self.__type = type
79
80 def get_min(self, v, o):
81 """
82 Parameters
83 ----------
84 v : initial value
85 o : object being fitted
86
87 Returns
88 -------
89 The minimum acceptable value for the range
90 """
91 return self.__callable(v, o)[0]
92
93 def get_max(self, v, o):
94 """
95 Parameters
96 ----------
97 v : initial value
98 o : object being fitted
99
100 Returns
101 -------
102 The maximum acceptable value for the range
103 """
104 return self.__callable(v, o)[1]
105
106 def get_type(self):
107 """
108 Returns
109 -------
110 RangeType
111 """
112 return self.__type
113
114 def __str__(self):
115 """
116 Returns
117 -------
118 str
119 Human readable representation for the object
120 """
121 res = '['
122 if hasattr(self.__limits, '__call__'):
123 res += 'func'
124 else:
125 res += '{},{}'.format(self.__limits[0], self.__limits[1])
126 type_str = {
127 RangeType.LINEAR: 'LIN',
128 RangeType.EXPONENTIAL: 'EXP'
129 }
130 res += ',{}]'.format(type_str[self.__type])
131 return res
132
133
134class Unbounded(object):
135 """
136 Unbounded, but normalize, value of a model fitting parameter
137
138 Parameters
139 ----------
140 normalization_factor: float, or a callable that receives the initial value parameter value and a source,
141 and returns a float
142 The world value which will be normalized to 1 in engine coordinates
143 """
144
145 def __init__(self, normalization_factor=1):
146 """
147 Constructor.
148 """
149 self.__normalization_factor = normalization_factor
150 if hasattr(normalization_factor, '__call__'):
151 self.__callable = normalization_factor
152 else:
153 self.__callable = lambda v, o: normalization_factor
154
156 """
157 Parameters
158 ----------
159 v : initial value
160 o : object being fitted
161
162 Returns
163 -------
164 float
165 The world value which will be normalized to 1 in engine coordinates
166 """
167 return self.__callable(v, o)
168
169 def __str__(self):
170 """
171 Returns
172 -------
173 str
174 Human readable representation for the object
175 """
176 res = '['
177 if hasattr(self.__normalization_factor, '__call__'):
178 res += 'func'
179 else:
180 res += '{}'.format(self.__normalization_factor)
181 res += ']'
182 return res
183
184
185class ParameterBase(cpp.Id):
186 """
187 Base class for all model fitting parameter types.
188 Can not be used directly.
189 """
190
191 def __str__(self):
192 """
193 Returns
194 -------
195 str
196 Human readable representation for the object
197 """
198 return '(ID:{})'.format(self.id)
199
200
202 """
203 A parameter with a single value that remains constant. It will not be fitted.
204
205 Parameters
206 ----------
207 value : float, or callable that receives a source and returns a float
208 Value for this parameter
209 """
210
211 def __init__(self, value):
212 """
213 Constructor.
214 """
215 ParameterBase.__init__(self)
216 self.__value = value
217 self.__callable = value if hasattr(value, '__call__') else lambda o: value
218
219 def get_value(self, o):
220 """
221 Parameters
222 ----------
223 o : object being fitted
224
225 Returns
226 -------
227 float
228 Value of the constant parameter
229 """
230 return self.__callable(o)
231
232 def __str__(self):
233 """
234 Returns
235 -------
236 str
237 Human readable representation for the object
238 """
239 res = ParameterBase.__str__(self)[:-1] + ', value:'
240 if hasattr(self.__value, '__call__'):
241 res += 'func'
242 else:
243 res += str(self.__value)
244 return res + ')'
245
246
248 """
249 A parameter that will be fitted by the model fitting engine.
250
251 Parameters
252 ----------
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.
257
258 See Also
259 --------
260 Unbounded
261 Range
262
263 Examples
264 --------
265 >>> sersic = FreeParameter(2.0, Range((1.0, 7.0), RangeType.LINEAR))
266 """
267
268 def __init__(self, init_value, range=Unbounded()):
269 """
270 Constructor.
271 """
272 ParameterBase.__init__(self)
273 self.__init_value = init_value
274 self.__init_call = init_value if hasattr(init_value, '__call__') else lambda o: init_value
275 self.__range = range
276
277 def get_init_value(self, o):
278 """
279 Parameters
280 ----------
281 o : object being fitted
282
283 Returns
284 -------
285 float
286 Initial value for the free parameter
287 """
288 return self.__init_call(o)
289
290 def get_range(self):
291 """
292 Returns
293 -------
294 Unbounded or Range
295 """
296 return self.__range
297
298 def __str__(self):
299 """
300 Returns
301 -------
302 str
303 Human readable representation for the object
304 """
305 res = ParameterBase.__str__(self)[:-1] + ', init:'
306 if hasattr(self.__init_value, '__call__'):
307 res += 'func'
308 else:
309 res += str(self.__init_value)
310 if self.__range:
311 res += ', range:' + str(self.__range)
312 return res + ')'
313
314
316 """
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
319
320 Parameters
321 ----------
322 func : callable
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.
327
328 Examples
329 --------
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)
333 """
334
335 def __init__(self, func, *params):
336 """
337 Constructor.
338 """
339 ParameterBase.__init__(self)
340 self.func = func
341 self.params = list(params)
342
343
345 """
346 Convenience function for the position parameter X and Y.
347
348 Returns
349 -------
350 x : FreeParameter
351 X coordinate, starting at the X coordinate of the centroid and linearly limited to X +/- the object radius.
352 y : FreeParameter
353 Y coordinate, starting at the Y coordinate of the centroid and linearly limited to Y +/- the object radius.
354 Notes
355 -----
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.
358 """
359 param_range = Range(lambda v, o: (v - o.radius, v + o.radius), RangeType.LINEAR)
360 return (
361 FreeParameter(lambda o: o.centroid_x, param_range),
362 FreeParameter(lambda o: o.centroid_y, param_range)
363 )
364
365
367 """
368 Possible flux types to use as initial value for the flux parameter.
369 Right now, only isophotal is supported.
370 """
371 ISO = 1
372
373
374def get_flux_parameter(type=FluxParameterType.ISO, scale=1):
375 """
376 Convenience function for the flux parameter.
377
378 Parameters
379 ----------
380 type : int
381 One of the values defined in FluxParameterType
382 scale : float
383 Scaling of the initial flux. Defaults to 1.
384
385 Returns
386 -------
387 flux : FreeParameter
388 Flux parameter, starting at the flux defined by `type`, and limited to +/- 1e3 times the initial value.
389 """
390 attr_map = {
391 FluxParameterType.ISO: 'isophotal_flux'
392 }
393 return FreeParameter(lambda o: getattr(o, attr_map[type]) * scale,
394 Range(lambda v, o: (v * 1E-3, v * 1E3), RangeType.EXPONENTIAL))
395
396
397class Prior(cpp.Id):
398 """
399 Model a Gaussian prior on a given parameter.
400
401 Parameters
402 ----------
403 param : ParameterBase
404 Model fitting parameter
405 value : float or callable that receives a source and returns a float
406 Mean of the Gaussian
407 sigma : float or callable that receives a source and returns a float
408 Standard deviation of the Gaussian
409 """
410
411 def __init__(self, param, value, sigma):
412 """
413 Constructor.
414 """
415 cpp.Id.__init__(self)
416 self.param = param.id
417 self.value = value if hasattr(value, '__call__') else lambda o: value
418 self.sigma = sigma if hasattr(sigma, '__call__') else lambda o: sigma
419
420
421class ModelBase(cpp.Id):
422 """
423 Base class for all models.
424 """
425
426 def get_params(self):
427 return []
428
429
431 """
432 Base class for positioned models with a flux. It can not be used directly.
433
434 Parameters
435 ----------
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
441 Total flux
442 """
443
444 def __init__(self, x_coord, y_coord, flux):
445 """
446 Constructor.
447 """
448 ModelBase.__init__(self)
449 self.x_coord = x_coord if isinstance(x_coord, ParameterBase) else ConstantParameter(x_coord)
450 self.y_coord = y_coord if isinstance(y_coord, ParameterBase) else ConstantParameter(y_coord)
451 self.flux = flux if isinstance(flux, ParameterBase) else ConstantParameter(flux)
452
453 def get_params(self):
454 return [self.x_coord, self.y_coord, self.flux]
455
456
457
459 """
460 Models a source as a point, spread by the PSF.
461
462 Parameters
463 ----------
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
469 Total flux
470 """
471
472 def __init__(self, x_coord, y_coord, flux):
473 """
474 Constructor.
475 """
476 CoordinateModelBase.__init__(self, x_coord, y_coord, flux)
477
478 def to_string(self, show_params=False):
479 """
480 Return a human readable representation of the model.
481
482 Parameters
483 ----------
484 show_params: bool
485 If True, include information about the parameters.
486
487 Returns
488 -------
489 str
490 """
491 if show_params:
492 return 'PointSource[x_coord={}, y_coord={}, flux={}]'.format(
493 self.x_coord, self.y_coord, self.flux)
494 else:
495 return 'PointSource[x_coord={}, y_coord={}, flux={}]'.format(
496 self.x_coord.id, self.y_coord.id, self.flux.id)
497
498
500 """
501 A model that is constant through all the image.
502
503 Parameters
504 ----------
505 value: ParameterBase or float
506 Value to add to the value of all pixels from the model.
507 """
508
509 def __init__(self, value):
510 """
511 Constructor.
512 """
513 ModelBase.__init__(self)
514 self.value = value if isinstance(value, ParameterBase) else ConstantParameter(value)
515
516 def to_string(self, show_params=False):
517 """
518 Return a human readable representation of the model.
519
520 Parameters
521 ----------
522 show_params: bool
523 If True, include information about the parameters.
524
525 Returns
526 -------
527 str
528 """
529 if show_params:
530 return 'ConstantModel[value={}]'.format(self.value)
531 else:
532 return 'ConstantModel[value={}]'.format(self.value.id)
533
534 def get_params(self):
535 return [self.value]
536
537
539 """
540 Base class for the Sersic, Exponential and de Vaucouleurs models. It can not be used directly.
541
542 Parameters
543 ----------
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
549 Total flux
550 effective_radius : ParameterBase or float
551 Ellipse semi-major axis, in pixels on the detection image.
552 aspect_ratio : ParameterBase or float
553 Ellipse ratio.
554 angle : ParameterBase or float
555 Ellipse rotation, in radians
556 """
557
558 def __init__(self, x_coord, y_coord, flux, effective_radius, aspect_ratio, angle):
559 """
560 Constructor.
561 """
562 CoordinateModelBase.__init__(self, x_coord, y_coord, flux)
563 self.effective_radius = effective_radius if isinstance(effective_radius,
564 ParameterBase) else ConstantParameter(
565 effective_radius)
566 self.aspect_ratio = aspect_ratio if isinstance(aspect_ratio,
567 ParameterBase) else ConstantParameter(
568 aspect_ratio)
569 self.angle = angle if isinstance(angle, ParameterBase) else ConstantParameter(angle)
570
571 def get_params(self):
572 return CoordinateModelBase.get_params(self) + [self.effective_radius, self.aspect_ratio, self.angle]
573
575 """
576 Model a source with a Sersic profile.
577
578 Parameters
579 ----------
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
585 Total flux
586 effective_radius : ParameterBase or float
587 Ellipse semi-major axis, in pixels on the detection image.
588 aspect_ratio : ParameterBase or float
589 Ellipse ratio.
590 angle : ParameterBase or float
591 Ellipse rotation, in radians
592 n : ParameterBase or float
593 Sersic index
594 """
595
596 def __init__(self, x_coord, y_coord, flux, effective_radius, aspect_ratio, angle, n):
597 """
598 Constructor.
599 """
600 SersicModelBase.__init__(self, x_coord, y_coord, flux, effective_radius, aspect_ratio,
601 angle)
602 self.n = n if isinstance(n, ParameterBase) else ConstantParameter(n)
603
604 def to_string(self, show_params=False):
605 """
606 Return a human readable representation of the model.
607
608 Parameters
609 ----------
610 show_params: bool
611 If True, include information about the parameters.
612
613 Returns
614 -------
615 str
616 """
617 if show_params:
618 return 'Sersic[x_coord={}, y_coord={}, flux={}, effective_radius={}, aspect_ratio={}, angle={}, n={}]'.format(
619 self.x_coord, self.y_coord, self.flux, self.effective_radius, self.aspect_ratio,
620 self.angle, self.n)
621 else:
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)
625
626 def get_params(self):
627 return SersicModelBase.get_params(self) + [self.n]
628
629
631 """
632 Model a source with an exponential profile (Sersic model with an index of 1)
633
634 Parameters
635 ----------
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
641 Total flux
642 effective_radius : ParameterBase or float
643 Ellipse semi-major axis, in pixels on the detection image.
644 aspect_ratio : ParameterBase or float
645 Ellipse ratio.
646 angle : ParameterBase or float
647 Ellipse rotation, in radians
648 """
649
650 def __init__(self, x_coord, y_coord, flux, effective_radius, aspect_ratio, angle):
651 """
652 Constructor.
653 """
654 SersicModelBase.__init__(self, x_coord, y_coord, flux, effective_radius, aspect_ratio,
655 angle)
656
657 def to_string(self, show_params=False):
658 """
659 Return a human readable representation of the model.
660
661 Parameters
662 ----------
663 show_params: bool
664 If True, include information about the parameters.
665
666 Returns
667 -------
668 str
669 """
670 if show_params:
671 return 'Exponential[x_coord={}, y_coord={}, flux={}, effective_radius={}, aspect_ratio={}, angle={}]'.format(
672 self.x_coord, self.y_coord, self.flux, self.effective_radius, self.aspect_ratio,
673 self.angle)
674 else:
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)
678
679
681 """
682 Model a source with a De Vaucouleurs profile (Sersic model with an index of 4)
683
684 Parameters
685 ----------
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
691 Total flux
692 effective_radius : ParameterBase or float
693 Ellipse semi-major axis, in pixels on the detection image.
694 aspect_ratio : ParameterBase or float
695 Ellipse ratio.
696 angle : ParameterBase or float
697 Ellipse rotation, in radians
698 """
699
700 def __init__(self, x_coord, y_coord, flux, effective_radius, aspect_ratio, angle):
701 """
702 Constructor.
703 """
704 SersicModelBase.__init__(self, x_coord, y_coord, flux, effective_radius, aspect_ratio,
705 angle)
706
707 def to_string(self, show_params=False):
708 """
709 Return a human readable representation of the model.
710
711 Parameters
712 ----------
713 show_params: bool
714 If True, include information about the parameters.
715
716 Returns
717 -------
718 str
719 """
720 if show_params:
721 return 'DeVaucouleurs[x_coord={}, y_coord={}, flux={}, effective_radius={}, aspect_ratio={}, angle={}]'.format(
722 self.x_coord, self.y_coord, self.flux, self.effective_radius, self.aspect_ratio,
723 self.angle)
724 else:
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)
728
729
731 """
732 ComputeGraphModel model
733
734 Parameters
735 ----------
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
743 Total flux
744 params : Dictionary of String and ParameterBase or float
745 Dictionary of custom parameters for the ONNX model
746 """
747
748 def __init__(self, models, x_coord, y_coord, flux, params={}):
749 """
750 Constructor.
751 """
752 CoordinateModelBase.__init__(self, x_coord, y_coord, flux)
753
754 ratio_name = "_aspect_ratio"
755 angle_name = "_angle"
756 scale_name = "_scale"
757
758 for k in params.keys():
759 if not isinstance(params[k], ParameterBase):
760 params[k] = ConstantParameter(params[k])
761
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
765
766 self.aspect_ratio = aspect_ratio if isinstance(aspect_ratio,
767 ParameterBase) else ConstantParameter(
768 aspect_ratio)
769 self.angle = angle if isinstance(angle, ParameterBase) else ConstantParameter(angle)
770 self.scale = scale if isinstance(scale, ParameterBase) else ConstantParameter(scale)
771
772 params.pop(ratio_name, None)
773 params.pop(angle_name, None)
774 params.pop(scale_name, None)
775
776 self.params = params
777 self.models = models if isinstance(models, list) else [models]
778
779 def to_string(self, show_params=False):
780 """
781 Return a human readable representation of the model.
782
783 Parameters
784 ----------
785 show_params: bool
786 If True, include information about the parameters.
787
788 Returns
789 -------
790 str
791 """
792 if show_params:
793 return 'ComputeGraph[x_coord={}, y_coord={}, flux={}, effective_radius={}, aspect_ratio={}, angle={}]'.format(
794 self.x_coord, self.y_coord, self.flux, self.effective_radius, self.aspect_ratio,
795 self.angle)
796 else:
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)
800
801 def get_params(self):
802 return (CoordinateModelBase.get_params(self) + [self.scale, self.aspect_ratio, self.angle] +
803 list(self.params.values()))
804
805
807 """
808 Coordinates in right ascension and declination
809
810 Parameters
811 ----------
812 ra : float
813 Right ascension
814 dec : float
815 Declination
816 """
817
818 def __init__(self, ra, dec):
819 """
820 Constructor.
821 """
822 self.ra = ra
823 self.dec = dec
824
825
827 """
828 Transform an (X, Y) in pixel coordinates on the detection image to (RA, DEC) in world coordinates.
829 Parameters
830 ----------
831 x : float
832 y : float
833
834 Returns
835 -------
836 WorldCoordinate
837 """
838 # -1 as we go FITS -> internal
839 wc_alpha = pyston.image_to_world_alpha(x - 1, y - 1)
840 wc_delta = pyston.image_to_world_delta(x - 1, y - 1)
841 return WorldCoordinate(wc_alpha, wc_delta)
842
843
845 """
846 Transform an (X, Y) in pixel coordinates on the detection image to astropy SkyCoord.
847
848 Parameters
849 ----------
850 x : float
851 y : float
852
853 Returns
854 -------
855 SkyCoord
856 """
857 coord = pixel_to_world_coordinate(x, y)
858 sky_coord = SkyCoord(ra=coord.ra * u.degree, dec=coord.dec * u.degree)
859 return sky_coord
860
861
862def radius_to_wc_angle(x, y, rad):
863 """
864 Transform a radius in pixels on the detection image to a radius in sky coordinates.
865
866 Parameters
867 ----------
868 x : float
869 y : float
870 rad : float
871
872 Returns
873 -------
874 Radius in degrees
875 """
876 return (get_separation_angle(x, y, x + rad, y) + get_separation_angle(x, y, x, y + rad)) / 2.0
877
878
879def get_separation_angle(x1, y1, x2, y2):
880 """
881 Get the separation angle in sky coordinates for two points defined in pixels on the detection image.
882
883 Parameters
884 ----------
885 x1 : float
886 y1 : float
887 x2 : float
888 y2 : float
889
890 Returns
891 -------
892 Separation in degrees
893 """
894 c1 = get_sky_coord(x1, y1)
895 c2 = get_sky_coord(x2, y2)
896 return c1.separation(c2).degree
897
898
899def get_position_angle(x1, y1, x2, y2):
900 """
901 Get the position angle in sky coordinates for two points defined in pixels on the detection image.
902
903 Parameters
904 ----------
905 x1
906 y1
907 x2
908 y2
909
910 Returns
911 -------
912 Position angle in degrees, normalized to -/+ 90
913 """
914 c1 = get_sky_coord(x1, y1)
915 c2 = get_sky_coord(x2, y2)
916 angle = c1.position_angle(c2).degree
917
918 # return angle normalized to range: -90 <= angle < 90
919 return (angle + 90.0) % 180.0 - 90.0
920
921
923 """
924 Convenience function for generating two dependent parameter with world (alpha, delta) coordinates
925 from image (X, Y) coordinates.
926
927 Parameters
928 ----------
929 x : ParameterBase
930 y : ParameterBase
931
932 Returns
933 -------
934 ra : DependentParameter
935 dec : DependentParameter
936
937 See Also
938 --------
939 get_pos_parameters
940
941 Examples
942 --------
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)
947 """
948 ra = DependentParameter(lambda x, y: pixel_to_world_coordinate(x, y).ra, x, y)
949 dec = DependentParameter(lambda x, y: pixel_to_world_coordinate(x, y).dec, x, y)
950 return (ra, dec)
951
952
953def get_world_parameters(x, y, radius, angle, ratio):
954 """
955 Convenience function for generating five dependent parameters, in world coordinates, for the position
956 and shape of a model.
957
958 Parameters
959 ----------
960 x : ParameterBase
961 y : ParameterBase
962 radius : ParameterBase
963 angle : ParameterBase
964 ratio : ParameterBase
965
966 Returns
967 -------
968 ra : DependentParameter
969 Right ascension
970 dec : DependentParameter
971 Declination
972 rad : DependentParameter
973 Radius as degrees
974 angle : DependentParameter
975 Angle in degrees
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
979
980 Examples
981 --------
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)
990 """
991 ra = DependentParameter(lambda x, y: pixel_to_world_coordinate(x, y).ra, x, y)
992 dec = DependentParameter(lambda x, y: pixel_to_world_coordinate(x, y).dec, x, y)
993
994 def get_major_axis(x, y, radius, angle, ratio):
995 if ratio <= 1:
996 x2 = x + math.cos(angle) * radius
997 y2 = y + math.sin(angle) * radius
998 else:
999 x2 = x + math.sin(angle) * radius * ratio
1000 y2 = y + math.cos(angle) * radius * ratio
1001
1002 return x2, y2
1003
1004 def get_minor_axis(x, y, radius, angle, ratio):
1005 if ratio <= 1:
1006 x2 = x + math.sin(angle) * radius * ratio
1007 y2 = y + math.cos(angle) * radius * ratio
1008 else:
1009 x2 = x + math.cos(angle) * radius
1010 y2 = y + math.sin(angle) * radius
1011
1012 return x2, y2
1013
1014 def wc_rad_func(x, y, radius, angle, ratio):
1015 x2, y2 = get_major_axis(x, y, radius, angle, ratio)
1016 return get_separation_angle(x, y, x2, y2)
1017
1018 def wc_angle_func(x, y, radius, angle, ratio):
1019 x2, y2 = get_major_axis(x, y, radius, angle, ratio)
1020 return get_position_angle(x, y, x2, y2)
1021
1022 def wc_ratio_func(x, y, radius, angle, ratio):
1023 minor_x, minor_y = get_minor_axis(x, y, radius, angle, ratio)
1024 minor_angle = get_separation_angle(x, y, minor_x, minor_y)
1025
1026 major_x, major_y = get_major_axis(x, y, radius, angle, ratio)
1027 major_angle = get_separation_angle(x, y, major_x, major_y)
1028
1029 return minor_angle / major_angle
1030
1031 wc_rad = DependentParameter(wc_rad_func, x, y, radius, angle, ratio)
1032 wc_angle = DependentParameter(wc_angle_func, x, y, radius, angle, ratio)
1033 wc_ratio = DependentParameter(wc_ratio_func, x, y, radius, angle, ratio)
1034
1035 return ra, dec, wc_rad, wc_angle, wc_ratio
1036
1037class WindowType(Enum):
1038 RECTANGLE = 1
1039 SQUARE_MIN = 2
1040 SQUARE_MAX = 3
1041 SQUARE_AREA = 4
1042 DISK_MIN = 5
1043 DISK_MAX = 6
1044 DISK_AREA = 7
1045 ALIGNED_ELLIPSE = 8
1046 ROTATED_ELLIPSE = 9
1047
1049 def __init__(self):
1054 self.prior_dict = {}
1061 self.params_dict = {"max_iterations": 200, "modified_chi_squared_scale": 10, "engine": "",
1062 "use_iterative_fitting": True, "meta_iterations": 5,
1063 "deblend_factor": 0.95, "meta_iteration_stop": 0.0001,
1064 "window_type": WindowType.RECTANGLE, "ellipse_scale": 3.0
1065 }
1066
1067 def _set_model_to_frames(self, group, model):
1068 for x in group:
1069 if isinstance(x, tuple):
1070 self._set_model_to_frames(x[1], model)
1071 else:
1072 if x.id not in self.frame_models_dict:
1073 self.frame_models_dict[x.id] = []
1074 self.frame_models_dict[x.id].append(model.id)
1075
1076 def _is_param_known(self, param):
1077 return param.id in self.constant_parameter_dict or \
1078 param.id in self.free_parameter_dict or \
1079 param.id in self.dependent_parameter_dict
1080
1081 def _register_parameter(self, attr):
1082 if isinstance(attr, ConstantParameter):
1083 self.constant_parameter_dict[attr.id] = attr
1084 elif isinstance(attr, FreeParameter):
1085 self.free_parameter_dict[attr.id] = attr
1086 elif isinstance(attr, DependentParameter):
1087 self.dependent_parameter_dict[attr.id] = attr
1088 for param in attr.params:
1089 self._register_parameter(param)
1090
1091 def _populate_parameters(self, model):
1092 for param in model.get_params():
1093 self._register_parameter(param)
1094
1095 def _register_model(self, model):
1096 if isinstance(model, ConstantModel):
1097 self.constant_model_dict[model.id] = model
1098 elif isinstance(model, PointSourceModel):
1099 self.point_source_model_dict[model.id] = model
1100 elif isinstance(model, SersicModel):
1101 self.sersic_model_dict[model.id] = model
1102 elif isinstance(model, ExponentialModel):
1103 self.exponential_model_dict[model.id] = model
1104 elif isinstance(model, DeVaucouleursModel):
1105 self.de_vaucouleurs_model_dict[model.id] = model
1106 elif isinstance(model, ComputeGraphModel):
1107 self.onnx_model_dict[model.id] = model
1108 else:
1109 raise TypeError('Unknown model type {}'.format(type(model)))
1110
1111 def add_model(self, group, model):
1112 """
1113 Add a model to be fitted to the given group.
1114
1115 Parameters
1116 ----------
1117 group : MeasurementGroup
1118 model : ModelBase
1119 """
1120 if not isinstance(group, MeasurementGroup):
1121 raise TypeError(
1122 'Models can only be added on MeasurementGroup, got {}'.format(type(group)))
1123 if not hasattr(group, 'models'):
1124 group.models = []
1125 group.models.append(model)
1126 self._set_model_to_frames(group, model)
1127 self._populate_parameters(model)
1128 self._register_model(model)
1129
1130 def add_prior(self, param, value, sigma):
1131 """
1132 Add a prior to the given parameter.
1133
1134 Parameters
1135 ----------
1136 param : ParameterBase
1137 value : float or callable that receives a source and returns a float
1138 Mean of the Gaussian
1139 sigma : float or callable that receives a source and returns a float
1140 Standard deviation of the Gaussian
1141 """
1142 prior = Prior(param, value, sigma)
1143 self.prior_dict[prior.id] = prior
1144 if not self._is_param_known(param):
1145 self._register_parameter(param)
1146
1147 def print_parameters(self, file=sys.stderr):
1148 """
1149 Print a human-readable representation of the configured model fitting parameters.
1150
1151 Parameters
1152 ----------
1153 file : file object
1154 Where to print the representation. Defaults to sys.stderr
1155 """
1157 print('Constant parameters:', file=file)
1158 for n, param in self.constant_parameter_dict.items():
1159 print(' {}: {}'.format(n, param), file=file)
1160 if self.free_parameter_dict:
1161 print('Free parameters:', file=file)
1162 for n, param in self.free_parameter_dict.items():
1163 print(' {}: {}'.format(n, param), file=file)
1165 print('Dependent parameters:', file=file)
1166 for n, param in self.dependent_parameter_dict.items():
1167 print(' {}: {}'.format(n, param), file=file)
1168
1169 def set_max_iterations(self, iterations):
1170 """
1171 Parameters
1172 ----------
1173 iterations : int
1174 Max number of iterations for the model fitting.
1175 """
1176 self.params_dict["max_iterations"] = iterations
1177
1179 """
1180 Parameters
1181 ----------
1182 scale : float
1183 Sets u0, as used by the modified chi squared residual comparator, a function that reduces the effect of large
1184 deviations.
1185 Refer to the SourceXtractor++ documentation for a better explanation of how residuals are computed and how
1186 this value affects the model fitting.
1187 """
1188 self.params_dict["modified_chi_squared_scale"] = scale
1189
1190 def set_engine(self, engine):
1191 """
1192 Parameters
1193 ----------
1194 engine : str
1195 Minimization engine for the model fitting : levmar or gsl
1196 """
1197 self.params_dict["engine"] = engine
1198
1199 def use_iterative_fitting(self, use_iterative_fitting):
1200 """
1201 Parameters
1202 ----------
1203 use_iterative_fitting : boolean
1204 use iterative model fitting or legacy
1205 """
1206 self.params_dict["use_iterative_fitting"] = use_iterative_fitting
1207
1208 def set_meta_iterations(self, meta_iterations):
1209 """
1210 Parameters
1211 ----------
1212 meta_iterations : int
1213 number of meta iterations on the whole group (when using iterative model fitting)
1214 """
1215 self.params_dict["meta_iterations"] = meta_iterations
1216
1217 def set_deblend_factor(self, deblend_factor):
1218 """
1219 Parameters
1220 ----------
1221
1222 """
1223 self.params_dict["deblend_factor"] = deblend_factor
1224
1225 def set_meta_iteration_stop(self, meta_iteration_stop):
1226 """
1227 Parameters
1228 ----------
1229
1230 """
1231 self.params_dict["meta_iteration_stop"] = meta_iteration_stop
1232
1233 def set_window_type(self, window_type):
1234 """
1235 Parameters
1236 ----------
1237
1238 window_type : WindowType
1239 specify the type of model fitting window
1240
1241 """
1242 if type(window_type) != WindowType:
1243 raise TypeError("Window type must be a WindowType enum value!")
1244
1245 self.params_dict["window_type"] = window_type
1246
1247 def set_ellipse_scale(self, ellipse_scale):
1248 """
1249 Parameters
1250 ----------
1251
1252 ellipse_scale : double
1253 specify scaling of elliptic footprints
1254
1255 """
1256 self.params_dict["ellipse_scale"] = ellipse_scale
1257
1258
1259def print_model_fitting_info(group, show_params=False, prefix='', file=sys.stderr):
1260 """
1261 Print a human-readable representation of the configured models.
1262
1263 Parameters
1264 ----------
1265 group : MeasurementGroup
1266 Print the models for this group.
1267 show_params : bool
1268 If True, print also the parameters that belong to the model
1269 prefix : str
1270 Prefix each line with this string. Used internally for indentation.
1271 file : file object
1272 Where to print the representation. Defaults to sys.stderr
1273 """
1274 if hasattr(group, 'models') and group.models:
1275 print('{}Models:'.format(prefix), file=file)
1276 for m in group.models:
1277 print('{} {}'.format(prefix, m.to_string(show_params)), file=file)
1278 for x in group:
1279 if isinstance(x, tuple):
1280 print('{}{}:'.format(prefix, x[0]), file=file)
1281 print_model_fitting_info(x[1], show_params, prefix + ' ', file=file)
__init__(self, models, x_coord, y_coord, flux, params={})
__init__(self, x_coord, y_coord, flux, effective_radius, aspect_ratio, angle)
__init__(self, x_coord, y_coord, flux, effective_radius, aspect_ratio, angle)
__init__(self, init_value, range=Unbounded())
use_iterative_fitting(self, use_iterative_fitting)
__init__(self, x_coord, y_coord, flux, effective_radius, aspect_ratio, angle)
__init__(self, x_coord, y_coord, flux, effective_radius, aspect_ratio, angle, n)
print_model_fitting_info(group, show_params=False, prefix='', file=sys.stderr)
get_world_parameters(x, y, radius, angle, ratio)
get_flux_parameter(type=FluxParameterType.ISO, scale=1)