LSSTApplications  16.0-10-g4f78f78+16,16.0-10-gc1446dd+42,16.0-11-g09ed895+1,16.0-13-g7649090,16.0-14-g0a28612+1,16.0-14-g6c7ed55+16,16.0-15-ga29f190+1,16.0-16-g89065d4+14,16.0-16-gd8e3590+16,16.0-16-ge6a35c8+6,16.0-17-g7e0e4ff+10,16.0-17-ga3d2e9f,16.0-19-gb830ed4e+16,16.0-2-g0febb12+21,16.0-2-g9d5294e+61,16.0-2-ga8830df+5,16.0-24-gc1c7f52+9,16.0-24-gd8faad8+3,16.0-3-ge00e371+21,16.0-36-g07840cb1,16.0-4-g18f3627+5,16.0-4-g5f3a788+20,16.0-4-ga3eb747+10,16.0-4-gabf74b7+16,16.0-4-gade8416+9,16.0-4-gb13d127+5,16.0-5-g6a53317+21,16.0-5-gb3f8a4b+74,16.0-5-gef99c9f+12,16.0-6-g9321be7+4,16.0-6-gcbc7b31+22,16.0-6-gf49912c+16,16.0-63-gae20905ba,16.0-7-gd2eeba5+31,16.0-8-g21fd5fe+16,16.0-8-g3a9f023+12,16.0-8-g4734f7a,16.0-8-ge230541+3,16.0-9-g85d1a16+16,master-g72da5ad0b5+6,w.2018.47
LSSTDataManagementBasePackage
Public Member Functions | Protected Attributes | List of all members
lsst::jointcal::TanPix2RaDec Class Reference

the transformation that handles pix to sideral transfos (Gnomonic, possibly with polynomial distortions). More...

#include <Gtransfo.h>

Inheritance diagram for lsst::jointcal::TanPix2RaDec:
lsst::jointcal::BaseTanWcs lsst::jointcal::Gtransfo

Public Member Functions

 TanPix2RaDec (GtransfoLin const &pix2Tan, Point const &tangentPoint, const GtransfoPoly *corrections=nullptr)
 pix2Tan describes the transfo from pix to tangent plane (degrees). More...
 
GtransfoPoly getPix2TangentPlane () const
 the transformation from pixels to tangent plane (degrees) More...
 
virtual void pix2TP (double xPixel, double yPixel, double &xTangentPlane, double &yTangentPlane) const
 transforms from pixel space to tangent plane (degrees) More...
 
 TanPix2RaDec ()
 
TanPix2RaDec operator* (GtransfoLin const &right) const
 composition with GtransfoLin More...
 
std::unique_ptr< GtransfocomposeAndReduce (GtransfoLin const &right) const
 Return a reduced composition of newTransfo = this(right()), or nullptr if it cannot be reduced. More...
 
TanRaDec2Pix inverted () const
 approximate inverse : it ignores corrections; More...
 
std::unique_ptr< GtransforoughInverse (const Frame &region) const
 Overload the "generic routine" (available for all Gtransfo types. More...
 
std::unique_ptr< GtransfoinverseTransfo (const double precision, const Frame &region) const
 Inverse transfo: returns a TanRaDec2Pix if there are no corrections, or the iterative solver if there are. More...
 
std::unique_ptr< Gtransfoclone () const
 returns a copy (allocated by new) of the transformation. More...
 
void dump (std::ostream &stream) const
 dumps the transfo coefficients to stream. More...
 
double fit (StarMatchList const &starMatchList)
 Not implemented yet, because we do it otherwise. More...
 
void apply (const double xIn, const double yIn, double &xOut, double &yOut) const
 Transform pixels to ICRS RA, Dec in degrees. More...
 
void apply (Point const &in, Point &out) const
 applies the tranfo to in and writes into out. Is indeed virtual. More...
 
Point apply (Point const &in) const
 All these apply(..) shadow the virtual one in derived classes, unless one writes "using Gtransfo::apply". More...
 
Frame apply (Frame const &inputframe, bool inscribed) const
 Transform a bounding box, taking either the inscribed or circumscribed box. More...
 
Point getTangentPoint () const
 Get the sky origin (CRVAL in FITS WCS terminology) in degrees. More...
 
GtransfoLin getLinPart () const
 The Linear part (corresponding to CD's and CRPIX's) More...
 
const GtransfoPolygetCorr () const
 Get a non-owning pointer to the correction transform polynomial. More...
 
void setCorrections (std::unique_ptr< GtransfoPoly > corrections)
 Assign the correction polynomial (what it means is left to derived classes) More...
 
Point getCrPix () const
 Get the pixel origin of the WCS (CRPIX in FITS WCS terminology, but zero-based) More...
 
std::string __str__ ()
 
void transformStar (FatPoint &in) const
 allows to write MyTransfo(MyStar) More...
 
virtual double getJacobian (Point const &point) const
 returns the local jacobian. More...
 
virtual double getJacobian (const double x, const double y) const
 returns the local jacobian. More...
 
virtual std::unique_ptr< GtransfocomposeAndReduce (Gtransfo const &right) const
 Return a reduced composition of newTransfo = this(right()), or nullptr if it cannot be reduced. More...
 
virtual void computeDerivative (Point const &where, GtransfoLin &derivative, const double step=0.01) const
 Computes the local Derivative of a transfo, w.r.t. More...
 
virtual GtransfoLin linearApproximation (Point const &where, const double step=0.01) const
 linear (local) approximation. More...
 
virtual void transformPosAndErrors (const FatPoint &in, FatPoint &out) const
 
virtual void transformErrors (Point const &where, const double *vIn, double *vOut) const
 transform errors (represented as double[3] in order V(xx),V(yy),Cov(xy)) More...
 
void getParams (double *params) const
 params should be at least Npar() long More...
 
void offsetParams (Eigen::VectorXd const &delta)
 
virtual double paramRef (const int i) const
 
virtual double & paramRef (const int i)
 
virtual void paramDerivatives (Point const &where, double *dx, double *dy) const
 Derivative w.r.t parameters. More...
 
virtual int getNpar () const
 returns the number of parameters (to compute chi2's) More...
 
virtual std::shared_ptr< ast::MappingtoAstMap (jointcal::Frame const &domain) const
 Create an equivalent AST mapping for this transformation, including an analytic inverse if possible. More...
 
void write (const std::string &fileName) const
 
virtual void write (std::ostream &stream) const
 

Protected Attributes

GtransfoLin linPix2Tan
 
std::unique_ptr< GtransfoPolycorr
 
double ra0
 
double dec0
 
double cos0
 
double sin0
 

Detailed Description

the transformation that handles pix to sideral transfos (Gnomonic, possibly with polynomial distortions).

Definition at line 600 of file Gtransfo.h.

Constructor & Destructor Documentation

◆ TanPix2RaDec() [1/2]

lsst::jointcal::TanPix2RaDec::TanPix2RaDec ( GtransfoLin const &  pix2Tan,
Point const &  tangentPoint,
const GtransfoPoly corrections = nullptr 
)

pix2Tan describes the transfo from pix to tangent plane (degrees).

TangentPoint in degrees. Corrections are applied between Lin and deprojection parts (as in Swarp).

Definition at line 1417 of file Gtransfo.cc.

1419  : BaseTanWcs(pix2Tan, tangentPoint, corrections) {}
BaseTanWcs(GtransfoLin const &pix2Tan, Point const &tangentPoint, const GtransfoPoly *corrections=nullptr)
Definition: Gtransfo.cc:1320

◆ TanPix2RaDec() [2/2]

lsst::jointcal::TanPix2RaDec::TanPix2RaDec ( )

Definition at line 1422 of file Gtransfo.cc.

1422 : BaseTanWcs(GtransfoLin(), Point(0, 0), nullptr) {}
BaseTanWcs(GtransfoLin const &pix2Tan, Point const &tangentPoint, const GtransfoPoly *corrections=nullptr)
Definition: Gtransfo.cc:1320

Member Function Documentation

◆ __str__()

std::string lsst::jointcal::Gtransfo::__str__ ( )
inlineinherited

Definition at line 94 of file Gtransfo.h.

94  {
96  dump(s);
97  return s.str();
98  }
T str(T... args)
solver_t * s
virtual void dump(std::ostream &stream=std::cout) const =0
dumps the transfo coefficients to stream.

◆ apply() [1/4]

void lsst::jointcal::Gtransfo::apply ( Point const &  in,
Point out 
) const
inlineinherited

applies the tranfo to in and writes into out. Is indeed virtual.

Definition at line 71 of file Gtransfo.h.

71 { apply(in.x, in.y, out.x, out.y); }
virtual void apply(const double xIn, const double yIn, double &xOut, double &yOut) const =0

◆ apply() [2/4]

Point lsst::jointcal::Gtransfo::apply ( Point const &  in) const
inlineinherited

All these apply(..) shadow the virtual one in derived classes, unless one writes "using Gtransfo::apply".

Definition at line 75 of file Gtransfo.h.

75  {
76  double xout, yout;
77  apply(in.x, in.y, xout, yout);
78  return Point(xout, yout);
79  }
virtual void apply(const double xIn, const double yIn, double &xOut, double &yOut) const =0

◆ apply() [3/4]

Frame lsst::jointcal::Gtransfo::apply ( Frame const &  inputframe,
bool  inscribed 
) const
inherited

Transform a bounding box, taking either the inscribed or circumscribed box.

Parameters
[in]inputframeThe frame to be transformed.
[in]inscribedReturn the inscribed (true) or circumscribed (false) box.
Returns
The transformed frame.

Definition at line 74 of file Gtransfo.cc.

74  {
75  // 2 opposite corners
76  double xtmin1, xtmax1, ytmin1, ytmax1;
77  apply(inputframe.xMin, inputframe.yMin, xtmin1, ytmin1);
78  apply(inputframe.xMax, inputframe.yMax, xtmax1, ytmax1);
79  Frame fr1(std::min(xtmin1, xtmax1), std::min(ytmin1, ytmax1), std::max(xtmin1, xtmax1),
80  std::max(ytmin1, ytmax1));
81  // 2 other corners
82  double xtmin2, xtmax2, ytmin2, ytmax2;
83  apply(inputframe.xMin, inputframe.yMax, xtmin2, ytmax2);
84  apply(inputframe.xMax, inputframe.yMin, xtmax2, ytmin2);
85  Frame fr2(std::min(xtmin2, xtmax2), std::min(ytmin2, ytmax2), std::max(xtmin2, xtmax2),
86  std::max(ytmin2, ytmax2));
87 
88  if (inscribed) return fr1 * fr2;
89  return fr1 + fr2;
90 }
T min(T... args)
T max(T... args)
virtual void apply(const double xIn, const double yIn, double &xOut, double &yOut) const =0

◆ apply() [4/4]

void lsst::jointcal::BaseTanWcs::apply ( const double  xIn,
const double  yIn,
double &  xOut,
double &  yOut 
) const
virtualinherited

Transform pixels to ICRS RA, Dec in degrees.

Implements lsst::jointcal::Gtransfo.

Definition at line 1351 of file Gtransfo.cc.

1351  {
1352  double l, m; // radians in the tangent plane
1353  pix2TP(xIn, yIn, l, m); // l, m in degrees.
1354  l = deg2rad(l);
1355  m = deg2rad(m); // now in radians
1356  // Code inspired from worldpos.c in wcssubs (ancestor of the wcslib)
1357  /* At variance with wcslib, it collapses the projection to a plane
1358  and expression of sidereal cooordinates into a single set of
1359  operations. */
1360  double dect = cos0 - m * sin0;
1361  if (dect == 0) {
1362  LOGL_WARN(_log, "No sidereal coordinates at pole!");
1363  xOut = 0;
1364  yOut = 0;
1365  return;
1366  }
1367  double rat = ra0 + atan2(l, dect);
1368  dect = atan(std::cos(rat - ra0) * (m * cos0 + sin0) / dect);
1369  if (rat - ra0 > M_PI) rat -= (2. * M_PI);
1370  if (rat - ra0 < -M_PI) rat += (2. * M_PI);
1371  if (rat < 0.0) rat += (2. * M_PI);
1372  // convert to degree
1373  xOut = rad2deg(rat);
1374  yOut = rad2deg(dect);
1375 }
T atan(T... args)
#define M_PI
Definition: ListMatch.cc:31
T atan2(T... args)
T cos(T... args)
#define LOGL_WARN(logger, message...)
Log a warn-level message using a varargs/printf style interface.
Definition: Log.h:545
int m
Definition: SpanSet.cc:49
virtual void pix2TP(double xPixel, double yPixel, double &xTangentPlane, double &yTangentPlane) const =0
Transform from pixels to tangent plane (degrees)

◆ clone()

std::unique_ptr< Gtransfo > lsst::jointcal::TanPix2RaDec::clone ( ) const
virtual

returns a copy (allocated by new) of the transformation.

Implements lsst::jointcal::Gtransfo.

Definition at line 1474 of file Gtransfo.cc.

1474  {
1476 }
GtransfoLin getLinPart() const
The Linear part (corresponding to CD&#39;s and CRPIX&#39;s)
Definition: Gtransfo.cc:1379
STL class.
Point getTangentPoint() const
Get the sky origin (CRVAL in FITS WCS terminology) in degrees.
Definition: Gtransfo.cc:1377
std::unique_ptr< GtransfoPoly > corr
Definition: Gtransfo.h:592

◆ composeAndReduce() [1/2]

std::unique_ptr< Gtransfo > lsst::jointcal::Gtransfo::composeAndReduce ( Gtransfo const &  right) const
virtualinherited

Return a reduced composition of newTransfo = this(right()), or nullptr if it cannot be reduced.

"Reduced" in this context means that they are capable of being merged into a single transform, for example, for two polynomials:

\[ f(x) = 1 + x^2, g(x) = -1 + 3x \]

we would have h = f.composeAndReduce(g) == 2 - 6x + 9x^2.

To be overloaded by derived classes if they can properly reduce the composition.

Parameters
rightThe transform to apply first.
Returns
The new reduced and composed gtransfo, or nullptr if no such reduction is possible.

Reimplemented in lsst::jointcal::GtransfoIdentity.

Definition at line 92 of file Gtransfo.cc.

93  { // by default no way to compose
94  return std::unique_ptr<Gtransfo>(nullptr);
95 }
STL class.

◆ composeAndReduce() [2/2]

std::unique_ptr< Gtransfo > lsst::jointcal::TanPix2RaDec::composeAndReduce ( GtransfoLin const &  right) const

Return a reduced composition of newTransfo = this(right()), or nullptr if it cannot be reduced.

"Reduced" in this context means that they are capable of being merged into a single transform, for example, for two polynomials:

\[ f(x) = 1 + x^2, g(x) = -1 + 3x \]

we would have h = f.composeAndReduce(g) == 2 - 6x + 9x^2.

To be overloaded by derived classes if they can properly reduce the composition.

Parameters
rightThe transform to apply first.
Returns
The new reduced and composed gtransfo, or nullptr if no such reduction is possible.

Definition at line 1424 of file Gtransfo.cc.

1424  {
1425  if (right.getOrder() == 1) {
1426  return std::make_unique<TanPix2RaDec>((*this) * (right));
1427  } else {
1428  return std::unique_ptr<Gtransfo>(nullptr);
1429  }
1430 }
T right(T... args)
STL class.

◆ computeDerivative()

void lsst::jointcal::Gtransfo::computeDerivative ( Point const &  where,
GtransfoLin derivative,
const double  step = 0.01 
) const
virtualinherited

Computes the local Derivative of a transfo, w.r.t.

the Derivative is represented by a GtransfoLin, in which (hopefully), the offset terms are zero.

position.

Step is used for numerical derivation.

Derivative should transform a vector of offsets into a vector of offsets.

Reimplemented in lsst::jointcal::GtransfoLin, lsst::jointcal::GtransfoPoly, and lsst::jointcal::GtransfoIdentity.

Definition at line 116 of file Gtransfo.cc.

116  {
117  double x = where.x;
118  double y = where.y;
119  double xp0, yp0;
120  apply(x, y, xp0, yp0);
121 
122  double xp, yp;
123  apply(x + step, y, xp, yp);
124  derivative.a11() = (xp - xp0) / step;
125  derivative.a21() = (yp - yp0) / step;
126  apply(x, y + step, xp, yp);
127  derivative.a12() = (xp - xp0) / step;
128  derivative.a22() = (yp - yp0) / step;
129  derivative.dx() = 0;
130  derivative.dy() = 0;
131 }
int y
Definition: SpanSet.cc:49
int const step
double x
virtual void apply(const double xIn, const double yIn, double &xOut, double &yOut) const =0

◆ dump()

void lsst::jointcal::TanPix2RaDec::dump ( std::ostream stream) const
virtual

dumps the transfo coefficients to stream.

Implements lsst::jointcal::Gtransfo.

Definition at line 1478 of file Gtransfo.cc.

1478  {
1479  stream << " TanPix2RaDec, lin part :" << endl << linPix2Tan;
1480  Point tp = getTangentPoint();
1481  stream << " tangent point " << tp.x << ' ' << tp.y << endl;
1482  Point crpix = getCrPix();
1483  stream << " crpix " << crpix.x << ' ' << crpix.y << endl;
1484  if (corr) stream << "PV correction: " << endl << *corr;
1485 }
table::PointKey< double > crpix
Definition: OldWcs.cc:131
T endl(T... args)
Point getCrPix() const
Get the pixel origin of the WCS (CRPIX in FITS WCS terminology, but zero-based)
Definition: Gtransfo.cc:1383
Point getTangentPoint() const
Get the sky origin (CRVAL in FITS WCS terminology) in degrees.
Definition: Gtransfo.cc:1377
std::unique_ptr< GtransfoPoly > corr
Definition: Gtransfo.h:592

◆ fit()

double lsst::jointcal::TanPix2RaDec::fit ( StarMatchList const &  starMatchList)
virtual

Not implemented yet, because we do it otherwise.

Implements lsst::jointcal::Gtransfo.

Definition at line 1487 of file Gtransfo.cc.

1487  {
1488  /* OK we could implement this routine, but it is
1489  probably useless since to do the match, we have to
1490  project from sky to tangent plane. When a match is
1491  found, it is easier to carry out the fit in the
1492  tangent plane, rather than going back to the celestial
1493  sphere (and reproject to fit...). Anyway if this
1494  message shows up, we'll think about it.
1495  */
1496  throw LSST_EXCEPT(pex::exceptions::InvalidParameterError,
1497  "TanPix2RaDec::fit is NOT implemented (although it is doable)) ");
1498  return -1;
1499 }
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition: Exception.h:48

◆ getCorr()

const GtransfoPoly* lsst::jointcal::BaseTanWcs::getCorr ( ) const
inlineinherited

Get a non-owning pointer to the correction transform polynomial.

Definition at line 572 of file Gtransfo.h.

572 { return corr.get(); }
std::unique_ptr< GtransfoPoly > corr
Definition: Gtransfo.h:592

◆ getCrPix()

Point lsst::jointcal::BaseTanWcs::getCrPix ( ) const
inherited

Get the pixel origin of the WCS (CRPIX in FITS WCS terminology, but zero-based)

Definition at line 1383 of file Gtransfo.cc.

1383  {
1384  /* CRPIX's are defined by:
1385  ( CD1_1 CD1_2 ) (x - crpix1)
1386  transformed = ( ) * ( )
1387  ( CD2_1 CD2_2 ) (y - crpix2)
1388 
1389  so that CrPix is the point which transforms to (0,0)
1390  */
1391  const GtransfoLin inverse = linPix2Tan.inverted();
1392  return Point(inverse.Dx(), inverse.Dy());
1393 }
GtransfoLin inverted() const
returns the inverse: T1 = T2.inverted();
Definition: Gtransfo.cc:1202

◆ getJacobian() [1/2]

virtual double lsst::jointcal::Gtransfo::getJacobian ( Point const &  point) const
inlinevirtualinherited

returns the local jacobian.

Definition at line 111 of file Gtransfo.h.

111 { return getJacobian(point.x, point.y); }
virtual double getJacobian(Point const &point) const
returns the local jacobian.
Definition: Gtransfo.h:111

◆ getJacobian() [2/2]

double lsst::jointcal::Gtransfo::getJacobian ( const double  x,
const double  y 
) const
virtualinherited

returns the local jacobian.

Definition at line 97 of file Gtransfo.cc.

97  {
98  double x2, y2;
99  double eps = x * 0.01;
100  if (eps == 0) eps = 0.01;
101  apply(x, y, x2, y2);
102  double dxdx, dydx;
103  apply(x + eps, y, dxdx, dydx);
104  dxdx -= x2;
105  dydx -= y2;
106  double dxdy, dydy;
107  apply(x, y + eps, dxdy, dydy);
108  dxdy -= x2;
109  dydy -= y2;
110  return ((dxdx * dydy - dxdy * dydx) / (eps * eps));
111 }
int y
Definition: SpanSet.cc:49
double x
virtual void apply(const double xIn, const double yIn, double &xOut, double &yOut) const =0

◆ getLinPart()

GtransfoLin lsst::jointcal::BaseTanWcs::getLinPart ( ) const
inherited

The Linear part (corresponding to CD's and CRPIX's)

Definition at line 1379 of file Gtransfo.cc.

1379 { return linPix2Tan; }

◆ getNpar()

virtual int lsst::jointcal::Gtransfo::getNpar ( ) const
inlinevirtualinherited

returns the number of parameters (to compute chi2's)

Reimplemented in lsst::jointcal::GtransfoLinScale, lsst::jointcal::GtransfoLinRot, lsst::jointcal::GtransfoLinShift, lsst::jointcal::GtransfoPoly, and lsst::jointcal::GtransfoIdentity.

Definition at line 180 of file Gtransfo.h.

180 { return 0; }

◆ getParams()

void lsst::jointcal::Gtransfo::getParams ( double *  params) const
inherited

params should be at least Npar() long

Definition at line 210 of file Gtransfo.cc.

210  {
211  int npar = getNpar();
212  for (int i = 0; i < npar; ++i) params[i] = paramRef(i);
213 }
virtual int getNpar() const
returns the number of parameters (to compute chi2&#39;s)
Definition: Gtransfo.h:180
virtual double paramRef(const int i) const
Definition: Gtransfo.cc:220

◆ getPix2TangentPlane()

GtransfoPoly lsst::jointcal::TanPix2RaDec::getPix2TangentPlane ( ) const
virtual

the transformation from pixels to tangent plane (degrees)

Implements lsst::jointcal::BaseTanWcs.

Definition at line 1457 of file Gtransfo.cc.

1457  {
1458  if (corr)
1459  return (*corr) * linPix2Tan;
1460  else
1461  return linPix2Tan;
1462 }
std::unique_ptr< GtransfoPoly > corr
Definition: Gtransfo.h:592

◆ getTangentPoint()

Point lsst::jointcal::BaseTanWcs::getTangentPoint ( ) const
inherited

Get the sky origin (CRVAL in FITS WCS terminology) in degrees.

Definition at line 1377 of file Gtransfo.cc.

1377 { return Point(rad2deg(ra0), rad2deg(dec0)); }

◆ inverseTransfo()

std::unique_ptr< Gtransfo > lsst::jointcal::TanPix2RaDec::inverseTransfo ( const double  precision,
const Frame region 
) const
virtual

Inverse transfo: returns a TanRaDec2Pix if there are no corrections, or the iterative solver if there are.

Reimplemented from lsst::jointcal::Gtransfo.

Definition at line 1450 of file Gtransfo.cc.

1450  {
1451  if (!corr)
1452  return std::unique_ptr<Gtransfo>(new TanRaDec2Pix(getLinPart().inverted(), getTangentPoint()));
1453  else
1454  return std::unique_ptr<Gtransfo>(new GtransfoInverse(this, precision, region));
1455 }
TanRaDec2Pix inverted() const
approximate inverse : it ignores corrections;
Definition: Gtransfo.cc:1438
GtransfoLin getLinPart() const
The Linear part (corresponding to CD&#39;s and CRPIX&#39;s)
Definition: Gtransfo.cc:1379
STL class.
Point getTangentPoint() const
Get the sky origin (CRVAL in FITS WCS terminology) in degrees.
Definition: Gtransfo.cc:1377
std::unique_ptr< GtransfoPoly > corr
Definition: Gtransfo.h:592

◆ inverted()

TanRaDec2Pix lsst::jointcal::TanPix2RaDec::inverted ( ) const

approximate inverse : it ignores corrections;

Definition at line 1438 of file Gtransfo.cc.

1438  {
1439  if (corr != nullptr) {
1440  LOGL_WARN(_log, "You are inverting a TanPix2RaDec with corrections.");
1441  LOGL_WARN(_log, "The inverse you get ignores the corrections!");
1442  }
1443  return TanRaDec2Pix(getLinPart().inverted(), getTangentPoint());
1444 }
TanRaDec2Pix inverted() const
approximate inverse : it ignores corrections;
Definition: Gtransfo.cc:1438
GtransfoLin getLinPart() const
The Linear part (corresponding to CD&#39;s and CRPIX&#39;s)
Definition: Gtransfo.cc:1379
#define LOGL_WARN(logger, message...)
Log a warn-level message using a varargs/printf style interface.
Definition: Log.h:545
Point getTangentPoint() const
Get the sky origin (CRVAL in FITS WCS terminology) in degrees.
Definition: Gtransfo.cc:1377
std::unique_ptr< GtransfoPoly > corr
Definition: Gtransfo.h:592

◆ linearApproximation()

GtransfoLin lsst::jointcal::Gtransfo::linearApproximation ( Point const &  where,
const double  step = 0.01 
) const
virtualinherited

linear (local) approximation.

Reimplemented in lsst::jointcal::GtransfoLin, and lsst::jointcal::GtransfoIdentity.

Definition at line 133 of file Gtransfo.cc.

133  {
134  Point outwhere = apply(where);
135  GtransfoLin der;
136  computeDerivative(where, der, step);
137  return GtransfoLinShift(outwhere.x, outwhere.y) * der * GtransfoLinShift(-where.x, -where.y);
138 }
virtual void computeDerivative(Point const &where, GtransfoLin &derivative, const double step=0.01) const
Computes the local Derivative of a transfo, w.r.t.
Definition: Gtransfo.cc:116
int const step
virtual void apply(const double xIn, const double yIn, double &xOut, double &yOut) const =0

◆ offsetParams()

void lsst::jointcal::Gtransfo::offsetParams ( Eigen::VectorXd const &  delta)
inherited

Definition at line 215 of file Gtransfo.cc.

215  {
216  int npar = getNpar();
217  for (int i = 0; i < npar; ++i) paramRef(i) += delta[i];
218 }
virtual int getNpar() const
returns the number of parameters (to compute chi2&#39;s)
Definition: Gtransfo.h:180
virtual double paramRef(const int i) const
Definition: Gtransfo.cc:220

◆ operator*()

TanPix2RaDec lsst::jointcal::TanPix2RaDec::operator* ( GtransfoLin const &  right) const

composition with GtransfoLin

Definition at line 1432 of file Gtransfo.cc.

1432  {
1433  TanPix2RaDec result(*this);
1434  result.linPix2Tan = result.linPix2Tan * right;
1435  return result;
1436 }
py::object result
Definition: schema.cc:284
T right(T... args)

◆ paramDerivatives()

void lsst::jointcal::Gtransfo::paramDerivatives ( Point const &  where,
double *  dx,
double *  dy 
) const
virtualinherited

Derivative w.r.t parameters.

Derivatives should be al least 2*NPar long. first Npar, for x, last Npar for y.

Reimplemented in lsst::jointcal::GtransfoPoly.

Definition at line 229 of file Gtransfo.cc.

229  {
230  throw LSST_EXCEPT(pex::exceptions::InvalidParameterError,
231  "Gtransfo::paramDerivatives() should never be called ");
232 }
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition: Exception.h:48

◆ paramRef() [1/2]

double lsst::jointcal::Gtransfo::paramRef ( const int  i) const
virtualinherited

Reimplemented in lsst::jointcal::GtransfoPoly.

Definition at line 220 of file Gtransfo.cc.

220  {
221  throw LSST_EXCEPT(pex::exceptions::InvalidParameterError,
222  std::string("Gtransfo::paramRef should never be called "));
223 }
STL class.
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition: Exception.h:48

◆ paramRef() [2/2]

double & lsst::jointcal::Gtransfo::paramRef ( const int  i)
virtualinherited

Reimplemented in lsst::jointcal::GtransfoPoly.

Definition at line 225 of file Gtransfo.cc.

225  {
226  throw LSST_EXCEPT(pex::exceptions::InvalidParameterError, "Gtransfo::paramRef should never be called ");
227 }
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition: Exception.h:48

◆ pix2TP()

void lsst::jointcal::TanPix2RaDec::pix2TP ( double  xPixel,
double  yPixel,
double &  xTangentPlane,
double &  yTangentPlane 
) const
virtual

transforms from pixel space to tangent plane (degrees)

Implements lsst::jointcal::BaseTanWcs.

Definition at line 1464 of file Gtransfo.cc.

1464  {
1465  // xTangentPlane, yTangentPlane in degrees.
1466  linPix2Tan.apply(xPixel, yPixel, xTangentPlane, yTangentPlane);
1467  if (corr) {
1468  double xtmp = xTangentPlane;
1469  double ytmp = yTangentPlane;
1470  corr->apply(xtmp, ytmp, xTangentPlane, yTangentPlane); // still in degrees.
1471  }
1472 }
std::unique_ptr< GtransfoPoly > corr
Definition: Gtransfo.h:592
void apply(const double xIn, const double yIn, double &xOut, double &yOut) const override
Definition: Gtransfo.cc:546

◆ roughInverse()

std::unique_ptr< Gtransfo > lsst::jointcal::TanPix2RaDec::roughInverse ( const Frame region) const
virtual

Overload the "generic routine" (available for all Gtransfo types.

Reimplemented from lsst::jointcal::Gtransfo.

Definition at line 1446 of file Gtransfo.cc.

1446  {
1447  return std::unique_ptr<Gtransfo>(new TanRaDec2Pix(getLinPart().inverted(), getTangentPoint()));
1448 }
TanRaDec2Pix inverted() const
approximate inverse : it ignores corrections;
Definition: Gtransfo.cc:1438
GtransfoLin getLinPart() const
The Linear part (corresponding to CD&#39;s and CRPIX&#39;s)
Definition: Gtransfo.cc:1379
STL class.
Point getTangentPoint() const
Get the sky origin (CRVAL in FITS WCS terminology) in degrees.
Definition: Gtransfo.cc:1377

◆ setCorrections()

void lsst::jointcal::BaseTanWcs::setCorrections ( std::unique_ptr< GtransfoPoly corrections)
inherited

Assign the correction polynomial (what it means is left to derived classes)

Definition at line 1381 of file Gtransfo.cc.

1381 { corr = std::move(corrections); }
T move(T... args)
std::unique_ptr< GtransfoPoly > corr
Definition: Gtransfo.h:592

◆ toAstMap()

virtual std::shared_ptr<ast::Mapping> lsst::jointcal::Gtransfo::toAstMap ( jointcal::Frame const &  domain) const
inlinevirtualinherited

Create an equivalent AST mapping for this transformation, including an analytic inverse if possible.

Parameters
domainThe domain of the transfo, to help find an inverse.
Returns
An AST Mapping that represents this transformation.

Reimplemented in lsst::jointcal::GtransfoPoly, and lsst::jointcal::GtransfoIdentity.

Definition at line 189 of file Gtransfo.h.

189  {
190  throw std::logic_error("toAstMap is not implemented for this class.");
191  }
STL class.

◆ transformErrors()

void lsst::jointcal::Gtransfo::transformErrors ( Point const &  where,
const double *  vIn,
double *  vOut 
) const
virtualinherited

transform errors (represented as double[3] in order V(xx),V(yy),Cov(xy))

Definition at line 157 of file Gtransfo.cc.

157  {
158  GtransfoLin der;
159  computeDerivative(where, der, 0.01);
160  double a11 = der.A11();
161  double a22 = der.A22();
162  double a21 = der.A21();
163  double a12 = der.A12();
164 
165  /* (a11 a12) (vxx vxy)
166  M = ( ) and V = ( )
167  (a21 a22) (xvy vyy)
168 
169  Vxx = Vin[0], vyy = Vin[1], Vxy = Vin[2];
170  we want to compute M*V*tp(M)
171  A lin alg light package would be perfect...
172  */
173  int xx = 0;
174  int yy = 1;
175  int xy = 2;
176  // M*V :
177 
178  double b11 = a11 * vIn[xx] + a12 * vIn[xy];
179  double b22 = a21 * vIn[xy] + a22 * vIn[yy];
180  double b12 = a11 * vIn[xy] + a12 * vIn[yy];
181  double b21 = a21 * vIn[xx] + a22 * vIn[xy];
182 
183  // (M*V) * tp(M)
184 
185  vOut[xx] = b11 * a11 + b12 * a12;
186  vOut[xy] = b11 * a21 + b12 * a22;
187  vOut[yy] = b21 * a21 + b22 * a22;
188 }
virtual void computeDerivative(Point const &where, GtransfoLin &derivative, const double step=0.01) const
Computes the local Derivative of a transfo, w.r.t.
Definition: Gtransfo.cc:116

◆ transformPosAndErrors()

void lsst::jointcal::Gtransfo::transformPosAndErrors ( const FatPoint in,
FatPoint out 
) const
virtualinherited

Reimplemented in lsst::jointcal::TanRaDec2Pix, and lsst::jointcal::GtransfoPoly.

Definition at line 140 of file Gtransfo.cc.

140  {
141  FatPoint res; // in case in and out are the same address...
142  res = apply(in);
143  GtransfoLin der;
144  // could save a call here, since Derivative needs the transform of where that we already have
145  // 0.01 may not be a very good idea in all cases. May be we should provide a way of altering that.
146  computeDerivative(in, der, 0.01);
147  double a11 = der.A11();
148  double a22 = der.A22();
149  double a21 = der.A21();
150  double a12 = der.A12();
151  res.vx = a11 * (a11 * in.vx + 2 * a12 * in.vxy) + a12 * a12 * in.vy;
152  res.vy = a21 * a21 * in.vx + a22 * a22 * in.vy + 2. * a21 * a22 * in.vxy;
153  res.vxy = a21 * a11 * in.vx + a22 * a12 * in.vy + (a21 * a12 + a11 * a22) * in.vxy;
154  out = res;
155 }
virtual void computeDerivative(Point const &where, GtransfoLin &derivative, const double step=0.01) const
Computes the local Derivative of a transfo, w.r.t.
Definition: Gtransfo.cc:116
virtual void apply(const double xIn, const double yIn, double &xOut, double &yOut) const =0

◆ transformStar()

void lsst::jointcal::Gtransfo::transformStar ( FatPoint in) const
inlineinherited

allows to write MyTransfo(MyStar)

Definition at line 108 of file Gtransfo.h.

108 { transformPosAndErrors(in, in); }
virtual void transformPosAndErrors(const FatPoint &in, FatPoint &out) const
Definition: Gtransfo.cc:140

◆ write() [1/2]

void lsst::jointcal::Gtransfo::write ( const std::string fileName) const
inherited

Definition at line 239 of file Gtransfo.cc.

239  {
240  ofstream s(fileName.c_str());
241  write(s);
242  bool ok = !s.fail();
243  s.close();
244  if (!ok)
245  throw LSST_EXCEPT(pex::exceptions::InvalidParameterError,
246  "Gtransfo::write, something went wrong for file " + fileName);
247 }
void write(const std::string &fileName) const
Definition: Gtransfo.cc:239
STL class.
solver_t * s
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition: Exception.h:48
T c_str(T... args)

◆ write() [2/2]

void lsst::jointcal::Gtransfo::write ( std::ostream stream) const
virtualinherited

Reimplemented in lsst::jointcal::GtransfoPoly, and lsst::jointcal::GtransfoIdentity.

Definition at line 249 of file Gtransfo.cc.

249  {
250  throw LSST_EXCEPT(pex::exceptions::InvalidParameterError,
251  "Gtransfo::write(ostream), should never be called. MEans that it is missing in some "
252  "derived class ");
253 }
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition: Exception.h:48

Member Data Documentation

◆ corr

std::unique_ptr<GtransfoPoly> lsst::jointcal::BaseTanWcs::corr
protectedinherited

Definition at line 592 of file Gtransfo.h.

◆ cos0

double lsst::jointcal::BaseTanWcs::cos0
protectedinherited

Definition at line 594 of file Gtransfo.h.

◆ dec0

double lsst::jointcal::BaseTanWcs::dec0
protectedinherited

Definition at line 593 of file Gtransfo.h.

◆ linPix2Tan

GtransfoLin lsst::jointcal::BaseTanWcs::linPix2Tan
protectedinherited

Definition at line 590 of file Gtransfo.h.

◆ ra0

double lsst::jointcal::BaseTanWcs::ra0
protectedinherited

Definition at line 593 of file Gtransfo.h.

◆ sin0

double lsst::jointcal::BaseTanWcs::sin0
protectedinherited

Definition at line 594 of file Gtransfo.h.


The documentation for this class was generated from the following files: