PCL
WorldTransformation.h
Go to the documentation of this file.
1 // ____ ______ __
2 // / __ \ / ____// /
3 // / /_/ // / / /
4 // / ____// /___ / /___ PixInsight Class Library
5 // /_/ \____//_____/ PCL 2.7.0
6 // ----------------------------------------------------------------------------
7 // pcl/WorldTransformation.h - Released 2024-06-18T15:48:54Z
8 // ----------------------------------------------------------------------------
9 // This file is part of the PixInsight Class Library (PCL).
10 // PCL is a multiplatform C++ framework for development of PixInsight modules.
11 //
12 // Copyright (c) 2003-2024 Pleiades Astrophoto S.L. All Rights Reserved.
13 //
14 // Redistribution and use in both source and binary forms, with or without
15 // modification, is permitted provided that the following conditions are met:
16 //
17 // 1. All redistributions of source code must retain the above copyright
18 // notice, this list of conditions and the following disclaimer.
19 //
20 // 2. All redistributions in binary form must reproduce the above copyright
21 // notice, this list of conditions and the following disclaimer in the
22 // documentation and/or other materials provided with the distribution.
23 //
24 // 3. Neither the names "PixInsight" and "Pleiades Astrophoto", nor the names
25 // of their contributors, may be used to endorse or promote products derived
26 // from this software without specific prior written permission. For written
27 // permission, please contact info@pixinsight.com.
28 //
29 // 4. All products derived from this software, in any form whatsoever, must
30 // reproduce the following acknowledgment in the end-user documentation
31 // and/or other materials provided with the product:
32 //
33 // "This product is based on software from the PixInsight project, developed
34 // by Pleiades Astrophoto and its contributors (https://pixinsight.com/)."
35 //
36 // Alternatively, if that is where third-party acknowledgments normally
37 // appear, this acknowledgment must be reproduced in the product itself.
38 //
39 // THIS SOFTWARE IS PROVIDED BY PLEIADES ASTROPHOTO AND ITS CONTRIBUTORS
40 // "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED
41 // TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
42 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL PLEIADES ASTROPHOTO OR ITS
43 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
44 // EXEMPLARY OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, BUSINESS
45 // INTERRUPTION; PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; AND LOSS OF USE,
46 // DATA OR PROFITS) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
47 // CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
48 // ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
49 // POSSIBILITY OF SUCH DAMAGE.
50 // ----------------------------------------------------------------------------
51 
52 #ifndef __PCL_WorldTransformation_h
53 #define __PCL_WorldTransformation_h
54 
56 
57 #include <pcl/Defs.h>
58 #include <pcl/Diagnostics.h>
59 
60 #include <pcl/GridInterpolation.h>
62 #include <pcl/SurfaceSpline.h>
63 #include <pcl/WCSKeywords.h>
64 
65 /*
66  * Based on original work contributed by AndrĂ©s del Pozo.
67  */
68 
69 namespace pcl
70 {
71 
72 // ----------------------------------------------------------------------------
73 
74 #define __PCL_WCS_DEFAULT_SPLINE_ORDER 2
75 #define __PCL_WCS_DEFAULT_SPLINE_SMOOTHNESS 0.005F
76 #define __PCL_WCS_SURFACE_SIMPLIFIER_DEFAULT_ENABLED true
77 #define __PCL_WCS_SURFACE_SIMPLIFIER_DEFAULT_REJECT_FRACTION 0.10F
78 #define __PCL_WCS_MAX_SPLINE_POINTS 2100
79 
80 // ----------------------------------------------------------------------------
81 
88 class PCL_CLASS WorldTransformation
89 {
90 public:
91 
95  WorldTransformation() = default;
96 
101 
106  {
107  }
108 
112  virtual bool IsEmpty() const
113  {
114  return false;
115  }
116 
120  virtual WorldTransformation* Clone() const = 0;
121 
129  virtual DPoint Direct( const DPoint& p ) const = 0;
130 
138  virtual DPoint Inverse( const DPoint& p ) const = 0;
139 
145 };
146 
147 // ----------------------------------------------------------------------------
148 
156 {
157 public:
158 
165  : m_transWI( transIW.Inverse() )
166  , m_transIW( transIW )
167  {
168  }
169 
174 
179 
184  {
185  }
186 
190  bool IsEmpty() const override
191  {
192  return false;
193  }
194 
197  WorldTransformation* Clone() const override
198  {
199  return new LinearWorldTransformation( *this );
200  }
201 
204  DPoint Direct( const DPoint& p ) const override
205  {
206  return m_transWI.Transform( p );
207  }
208 
211  DPoint Inverse( const DPoint& p ) const override
212  {
213  return m_transIW.Transform( p );
214  }
215 
221  {
222  return m_transIW;
223  }
224 
225 private:
226 
227  LinearTransformation m_transWI; // world -> image
228  LinearTransformation m_transIW; // image -> world
229 };
230 
231 // ----------------------------------------------------------------------------
232 
233 class PCL_CLASS XMLElement;
234 
251 {
252 public:
253 
340  SplineWorldTransformation( const Array<DPoint>& controlPointsW,
341  const Array<DPoint>& controlPointsI,
342  float smoothness = __PCL_WCS_DEFAULT_SPLINE_SMOOTHNESS,
343  const FVector& weights = FVector(),
344  int order = __PCL_WCS_DEFAULT_SPLINE_ORDER,
345  bool enableSimplifier = __PCL_WCS_SURFACE_SIMPLIFIER_DEFAULT_ENABLED,
346  float simplifierRejectFraction = __PCL_WCS_SURFACE_SIMPLIFIER_DEFAULT_REJECT_FRACTION )
347  : m_controlPointsW( controlPointsW )
348  , m_controlPointsI( controlPointsI )
349  , m_order( order )
350  , m_smoothness( smoothness )
351  , m_weights( weights )
352  , m_useSimplifiers( enableSimplifier )
353  , m_simplifierRejectFraction( simplifierRejectFraction )
354  {
355  InitializeSplines();
356  CalculateLinearApproximation();
357  }
358 
373  const LinearTransformation& linearIW = LinearTransformation::Null() );
374 
384  {
385  Deserialize( data );
386  InitializeSplines();
387  if ( linearIW.IsSingularMatrix() )
388  CalculateLinearApproximation();
389  else
390  m_linearIW = linearIW;
391  }
392 
397 
402 
407  {
408  }
409 
413  SplineWorldTransformation& operator =( const SplineWorldTransformation& ) = default;
414 
419 
426  bool IsEmpty() const override
427  {
428  return m_controlPointsW.IsEmpty() || m_controlPointsI.IsEmpty();
429  }
430 
434  WorldTransformation* Clone() const override
435  {
436  return new SplineWorldTransformation( *this );
437  }
438 
455  DPoint Direct( const DPoint& pW ) const override
456  {
457  if ( m_gridWI.IsValid() )
458  if ( m_gridWI.ReferenceRect().IncludesFast( pW ) )
459  return m_gridWI( pW );
460  return m_splineWI( pW );
461  }
462 
479  DPoint Inverse( const DPoint& pI ) const override
480  {
481  if ( m_gridIW.IsValid() )
482  if ( m_gridIW.ReferenceRect().IncludesFast( pI ) )
483  return m_gridIW( pI );
484  return m_splineIW( pI );
485  }
486 
492  {
493  return m_linearIW;
494  }
495 
500  const IsoString& Version() const
501  {
502  return m_version;
503  }
504 
521  void InitializeGridInterpolations( const Rect& rectI, int deltaI = 16 )
522  {
523  PCL_PRECONDITION( deltaI >= 1 && deltaI <= 64 )
524  deltaI = Range( deltaI, 1, 64 );
525 
526  /*
527  * Compute reference region and grid distance in native spherical
528  * coordinates.
529  */
530  DRect rectW( std::numeric_limits<double>::max(), std::numeric_limits<double>::max(),
531  std::numeric_limits<double>::lowest(), std::numeric_limits<double>::lowest() );
532  for ( int x = 0; x < rectI.Width(); ++x )
533  for ( int p = 0; p < 2; ++p )
534  {
535  DPoint pW = m_splineIW( DPoint( double( x ), p*(rectI.Height() - 1) ) );
536  if ( pW.x < rectW.x0 )
537  rectW.x0 = pW.x;
538  if ( pW.x > rectW.x1 )
539  rectW.x1 = pW.x;
540  if ( pW.y < rectW.y0 )
541  rectW.y0 = pW.y;
542  if ( pW.y > rectW.y1 )
543  rectW.y1 = pW.y;
544  }
545  for ( int y = 0; y < rectI.Height(); ++y )
546  for ( int p = 0; p < 2; ++p )
547  {
548  DPoint pW = m_splineIW( DPoint( p*(rectI.Width() - 1), double( y ) ) );
549  if ( pW.x < rectW.x0 )
550  rectW.x0 = pW.x;
551  if ( pW.x > rectW.x1 )
552  rectW.x1 = pW.x;
553  if ( pW.y < rectW.y0 )
554  rectW.y0 = pW.y;
555  if ( pW.y > rectW.y1 )
556  rectW.y1 = pW.y;
557  }
558  double deltaW = m_splineIW( DPoint( -deltaI/2, -deltaI/2 ) ).DistanceTo(
559  m_splineIW( DPoint( +deltaI/2, +deltaI/2 ) ) )/1.4142;
560 
561  m_gridWI.Initialize( rectW, deltaW, m_splineWI, false/*verbose*/ );
562  m_gridIW.Initialize( rectI, deltaI, m_splineIW, false/*verbose*/ );
563  }
564 
571  {
572  return m_gridWI.IsValid() && m_gridIW.IsValid();
573  }
574 
580  {
581  return "PCL:AstrometricSolution:SplineWorldTransformation:";
582  }
583 
589 
594  void Serialize( ByteArray& data ) const;
595 
608  {
609  return int( m_controlPointsW.Length() );
610  }
611 
625  {
626  return m_controlPointsW;
627  }
628 
642  {
643  return m_controlPointsI;
644  }
645 
650  int SplineOrder() const
651  {
652  return m_order;
653  }
654 
673  void GetSplineLengths( int& xWI, int& yWI, int& xIW, int& yIW ) const
674  {
675  xWI = m_splineWI.SplineX().Length();
676  yWI = m_splineWI.SplineY().Length();
677  xIW = m_splineIW.SplineX().Length();
678  yIW = m_splineIW.SplineY().Length();
679  }
680 
695  {
696  return m_truncated; // fetched during initialization
697  }
698 
704  bool SimplifiersEnabled() const
705  {
706  return m_useSimplifiers;
707  }
708 
715  {
716  return m_simplifierRejectFraction;
717  }
718 
734  void GetSplineErrors( double& xWI, double& yWI, double& xIW, double& yIW ) const
735  {
736  xWI = m_splineWI.ErrorX();
737  yWI = m_splineWI.ErrorY();
738  xIW = m_splineIW.ErrorX();
739  yIW = m_splineIW.ErrorY();
740  }
741 
742 private:
743 
744  IsoString m_version = "1.3";
745  Array<DPoint> m_controlPointsW;
746  Array<DPoint> m_controlPointsI;
747  int m_order = __PCL_WCS_DEFAULT_SPLINE_ORDER;
748  float m_smoothness = __PCL_WCS_DEFAULT_SPLINE_SMOOTHNESS;
749  FVector m_weights;
750  bool m_useSimplifiers = __PCL_WCS_SURFACE_SIMPLIFIER_DEFAULT_ENABLED;
751  float m_simplifierRejectFraction = __PCL_WCS_SURFACE_SIMPLIFIER_DEFAULT_REJECT_FRACTION;
752  bool m_truncated = false; // true => truncated control point vectors
753  Homography m_projectiveWI;
754  Homography m_projectiveIW;
755  PointSurfaceSpline m_splineWI;
756  PointSurfaceSpline m_splineIW;
757  PointGridInterpolation m_gridWI;
758  PointGridInterpolation m_gridIW;
759  LinearTransformation m_linearIW;
760 
761  static PointSurfaceSpline PointSurfaceSplineFromProperties( const PropertyArray& );
762  static PropertyArray PointSurfaceSplineToProperties( const PointSurfaceSpline&, const IsoString& );
763 
764  static SurfaceSpline<double> SurfaceSplineFromProperties( const PropertyArray& );
765  static PropertyArray SurfaceSplineToProperties( const SurfaceSpline<double>&, const IsoString& );
766 
767  static PointGridInterpolation PointGridInterpolationFromProperties( const PropertyArray& );
768  static PropertyArray PointGridInterpolationToProperties( const PointGridInterpolation&, const IsoString& );
769 
770  void Deserialize( const ByteArray& );
771  void InitializeSplines();
772  void CalculateLinearApproximation();
773 };
774 
775 // ----------------------------------------------------------------------------
776 
777 } // pcl
778 
779 #endif // __PCL_WorldTransformation_h
780 
781 // ----------------------------------------------------------------------------
782 // EOF pcl/WorldTransformation.h - Released 2024-06-18T15:48:54Z
A generic point in the two-dimensional space.
Definition: Point.h:100
component x
Abscissa (horizontal, or X-axis coordinate).
Definition: Point.h:111
component y
Ordinate (vertical, or Y-axis coordinate).
Definition: Point.h:112
A generic rectangle in the two-dimensional space.
Definition: Rectangle.h:314
component x1
Horizontal coordinate of the lower right corner.
Definition: Rectangle.h:334
component y1
Vertical coordinate of the lower right corner.
Definition: Rectangle.h:335
component y0
Vertical coordinate of the upper left corner.
Definition: Rectangle.h:333
component x0
Horizontal coordinate of the upper left corner.
Definition: Rectangle.h:332
component Width() const noexcept
Definition: Rectangle.h:635
component Height() const noexcept
Definition: Rectangle.h:644
Generic vector of arbitrary length.
Definition: Vector.h:107
Homography geometric transformation.
Definition: Homography.h:85
Eight-bit string (ISO/IEC-8859-1 or UTF-8 string)
Definition: String.h:5425
A linear geometric transformation on the plane defined as a 2x3 matrix of 64-bit floating point scala...
static LinearTransformation Null()
WCS linear world coordinate transformation.
WorldTransformation * Clone() const override
LinearWorldTransformation(LinearWorldTransformation &&)=default
LinearWorldTransformation(const LinearWorldTransformation &)=default
DPoint Direct(const DPoint &p) const override
DPoint Inverse(const DPoint &p) const override
LinearWorldTransformation(const LinearTransformation &transIW)
const LinearTransformation & ApproximateLinearTransform() const override
Discretized vector surface interpolation/approximation in two dimensions.
Vector surface spline interpolation/approximation in two dimensions.
Surface spline world coordinate transformation.
const Array< DPoint > & NativeControlPoints() const
DPoint Direct(const DPoint &pW) const override
void GetSplineErrors(double &xWI, double &yWI, double &xIW, double &yIW) const
SplineWorldTransformation(const ByteArray &data, const LinearTransformation &linearIW=LinearTransformation::Null())
WorldTransformation * Clone() const override
SplineWorldTransformation(const Array< DPoint > &controlPointsW, const Array< DPoint > &controlPointsI, float smoothness=__PCL_WCS_DEFAULT_SPLINE_SMOOTHNESS, const FVector &weights=FVector(), int order=__PCL_WCS_DEFAULT_SPLINE_ORDER, bool enableSimplifier=__PCL_WCS_SURFACE_SIMPLIFIER_DEFAULT_ENABLED, float simplifierRejectFraction=__PCL_WCS_SURFACE_SIMPLIFIER_DEFAULT_REJECT_FRACTION)
const IsoString & Version() const
const Array< DPoint > & ImageControlPoints() const
void Serialize(ByteArray &data) const
SplineWorldTransformation(const PropertyArray &properties, const LinearTransformation &linearIW=LinearTransformation::Null())
PropertyArray ToProperties() const
SplineWorldTransformation(const SplineWorldTransformation &)=default
void GetSplineLengths(int &xWI, int &yWI, int &xIW, int &yIW) const
void InitializeGridInterpolations(const Rect &rectI, int deltaI=16)
DPoint Inverse(const DPoint &pI) const override
SplineWorldTransformation(SplineWorldTransformation &&)=default
const LinearTransformation & ApproximateLinearTransform() const override
Abstract base class of world coordinate transformations.
WorldTransformation(const WorldTransformation &)=default
virtual bool IsEmpty() const
virtual WorldTransformation * Clone() const =0
virtual DPoint Inverse(const DPoint &p) const =0
virtual const LinearTransformation & ApproximateLinearTransform() const =0
virtual DPoint Direct(const DPoint &p) const =0
constexpr const T & Range(const T &x, const T &a, const T &b) noexcept
Definition: Utility.h:190
PCL root namespace.
Definition: AbstractImage.h:77