LSSTApplications  18.0.0+106,18.0.0+50,19.0.0,19.0.0+1,19.0.0+10,19.0.0+11,19.0.0+13,19.0.0+17,19.0.0+2,19.0.0-1-g20d9b18+6,19.0.0-1-g425ff20,19.0.0-1-g5549ca4,19.0.0-1-g580fafe+6,19.0.0-1-g6fe20d0+1,19.0.0-1-g7011481+9,19.0.0-1-g8c57eb9+6,19.0.0-1-gb5175dc+11,19.0.0-1-gdc0e4a7+9,19.0.0-1-ge272bc4+6,19.0.0-1-ge3aa853,19.0.0-10-g448f008b,19.0.0-12-g6990b2c,19.0.0-2-g0d9f9cd+11,19.0.0-2-g3d9e4fb2+11,19.0.0-2-g5037de4,19.0.0-2-gb96a1c4+3,19.0.0-2-gd955cfd+15,19.0.0-3-g2d13df8,19.0.0-3-g6f3c7dc,19.0.0-4-g725f80e+11,19.0.0-4-ga671dab3b+1,19.0.0-4-gad373c5+3,19.0.0-5-ga2acb9c+2,19.0.0-5-gfe96e6c+2,w.2020.01
LSSTDataManagementBasePackage
Public Member Functions | Public Attributes | List of all members
lsst::jointcal::Associations Class Reference

The class that implements the relations between MeasuredStar and FittedStar. More...

#include <Associations.h>

Public Member Functions

size_t refStarListSize ()
 
size_t fittedStarListSize ()
 
 Associations ()
 Source selection is performed in python, so Associations' constructor only initializes a couple of variables. More...
 
 Associations (CcdImageList const &imageList)
 Create an Associations object from a pre-built list of ccdImages. More...
 
 Associations (Associations const &)=delete
 No moves or copies: jointcal only ever needs one Associations object. More...
 
 Associations (Associations &&)=delete
 
Associationsoperator= (Associations const &)=delete
 
Associationsoperator= (Associations &&)=delete
 
void computeCommonTangentPoint ()
 Sets a shared tangent point for all ccdImages, using the mean of the centers of all ccdImages. More...
 
void setCommonTangentPoint (lsst::geom::Point2D const &commonTangentPoint)
 Sets a shared tangent point for all ccdImages. More...
 
Point getCommonTangentPoint () const
 can be used to project sidereal coordinates related to the image set on a plane. More...
 
void createCcdImage (afw::table::SourceCatalog &catalog, std::shared_ptr< lsst::afw::geom::SkyWcs > wcs, std::shared_ptr< lsst::afw::image::VisitInfo > visitInfo, lsst::geom::Box2I const &bbox, std::string const &filter, std::shared_ptr< afw::image::PhotoCalib > photoCalib, std::shared_ptr< afw::cameraGeom::Detector > detector, int visit, int ccd, lsst::jointcal::JointcalControl const &control)
 Create a ccdImage from an exposure catalog and metadata, and add it to the list. More...
 
void addCcdImage (std::shared_ptr< CcdImage > const ccdImage)
 Add a pre-constructed ccdImage to the ccdImageList. More...
 
void associateCatalogs (const double matchCutInArcsec=0, const bool useFittedList=false, const bool enlargeFittedList=true)
 incrementaly builds a merged catalog of all image catalogs More...
 
void collectRefStars (afw::table::SimpleCatalog &refCat, geom::Angle matchCut, std::string const &fluxField, float refCoordinateErr, bool rejectBadFluxes=false)
 Collect stars from an external reference catalog and associate them with fittedStars. More...
 
void deprojectFittedStars ()
 Sends back the fitted stars coordinates on the sky FittedStarsList::inTangentPlaneCoordinates keeps track of that. More...
 
void prepareFittedStars (int minMeasurements)
 Set the color field of FittedStar 's from a colored catalog. More...
 
CcdImageList const & getCcdImageList () const
 
unsigned getNFilters () const
 Number of different bands in the input image list. Not implemented so far. More...
 
lsst::sphgeom::Circle computeBoundingCircle () const
 Return the bounding circle in on-sky (RA, Dec) coordinates containing all CcdImages. More...
 
int nCcdImagesValidForFit () const
 return the number of CcdImages with non-empty catalogs to-be-fit. More...
 
size_t nFittedStarsWithAssociatedRefStar () const
 Return the number of fittedStars that have an associated refStar. More...
 

Public Attributes

CcdImageList ccdImageList
 
RefStarList refStarList
 
FittedStarList fittedStarList
 

Detailed Description

The class that implements the relations between MeasuredStar and FittedStar.

Definition at line 54 of file Associations.h.

Constructor & Destructor Documentation

◆ Associations() [1/4]

lsst::jointcal::Associations::Associations ( )
inline

Source selection is performed in python, so Associations' constructor only initializes a couple of variables.

Definition at line 69 of file Associations.h.

◆ Associations() [2/4]

lsst::jointcal::Associations::Associations ( CcdImageList const &  imageList)
inline

Create an Associations object from a pre-built list of ccdImages.

This is primarily useful for tests that build their own ccdImageList, but it could be used to help parallelize the creation of the ccdImages.

Parameters
imageListA pre-built ccdImage list.

Definition at line 81 of file Associations.h.

◆ Associations() [3/4]

lsst::jointcal::Associations::Associations ( Associations const &  )
delete

No moves or copies: jointcal only ever needs one Associations object.

◆ Associations() [4/4]

lsst::jointcal::Associations::Associations ( Associations &&  )
delete

Member Function Documentation

◆ addCcdImage()

void lsst::jointcal::Associations::addCcdImage ( std::shared_ptr< CcdImage > const  ccdImage)
inline

Add a pre-constructed ccdImage to the ccdImageList.

Definition at line 131 of file Associations.h.

131 { ccdImageList.push_back(ccdImage); }
T push_back(T... args)

◆ associateCatalogs()

void lsst::jointcal::Associations::associateCatalogs ( const double  matchCutInArcsec = 0,
const bool  useFittedList = false,
const bool  enlargeFittedList = true 
)

incrementaly builds a merged catalog of all image catalogs

Definition at line 121 of file Associations.cc.

122  {
123  // clear reference stars
124  refStarList.clear();
125 
126  // clear measurement counts and associations to refstars, but keep fittedStars themselves.
127  for (auto &item : fittedStarList) {
128  item->clearBeforeAssoc();
129  }
130  // clear fitted stars
131  if (!useFittedList) fittedStarList.clear();
132 
133  for (auto &ccdImage : ccdImageList) {
134  std::shared_ptr<AstrometryTransform> toCommonTangentPlane = ccdImage->getPixelToCommonTangentPlane();
135 
136  // Clear the catalog to fit and copy the whole catalog into it.
137  // This allows reassociating from scratch after a fit.
138  ccdImage->resetCatalogForFit();
139  MeasuredStarList &catalog = ccdImage->getCatalogForFit();
140 
141  // Associate with previous lists.
142  /* To speed up the match (more precisely the contruction of the FastFinder), select in the
143  fittedStarList the objects that are within reach of the current ccdImage */
144  Frame ccdImageFrameCPT = toCommonTangentPlane->apply(ccdImage->getImageFrame(), false);
145  ccdImageFrameCPT = ccdImageFrameCPT.rescale(1.10); // add 10 % margin.
146  // We cannot use FittedStarList::ExtractInFrame, because it does an actual copy, which we don't want
147  // here: we want the pointers in the StarMatch to refer to fittedStarList elements.
148  FittedStarList toMatch;
149 
150  for (auto const &fittedStar : fittedStarList) {
151  if (ccdImageFrameCPT.inFrame(*fittedStar)) {
152  toMatch.push_back(fittedStar);
153  }
154  }
155 
156  // divide by 3600 because coordinates in CTP are in degrees.
157  auto starMatchList = listMatchCollect(Measured2Base(catalog), Fitted2Base(toMatch),
158  toCommonTangentPlane.get(), matchCutInArcSec / 3600.);
159 
160  /* should check what this removeAmbiguities does... */
161  LOGLS_DEBUG(_log, "Measured-to-Fitted matches before removing ambiguities " << starMatchList->size());
162  starMatchList->removeAmbiguities(*toCommonTangentPlane);
163  LOGLS_DEBUG(_log, "Measured-to-Fitted matches after removing ambiguities " << starMatchList->size());
164 
165  // Associate MeasuredStar -> FittedStar using the surviving matches.
166 
167  int matchedCount = 0;
168  for (auto const &starMatch : *starMatchList) {
169  auto bs = starMatch.s1;
170  auto ms_const = std::dynamic_pointer_cast<const MeasuredStar>(bs);
171  auto ms = std::const_pointer_cast<MeasuredStar>(ms_const);
172  auto bs2 = starMatch.s2;
173  auto fs_const = std::dynamic_pointer_cast<const FittedStar>(bs2);
174  auto fs = std::const_pointer_cast<FittedStar>(fs_const);
175  ms->setFittedStar(fs);
176  matchedCount++;
177  }
178  LOGLS_INFO(_log, "Matched " << matchedCount << " objects in " << ccdImage->getName());
179 
180  // add unmatched objets to FittedStarList
181  int unMatchedCount = 0;
182  for (auto const &mstar : catalog) {
183  // to check if it was matched, just check if it has a fittedStar Pointer assigned
184  if (mstar->getFittedStar()) continue;
185  if (enlargeFittedList) {
186  auto fs = std::make_shared<FittedStar>(*mstar);
187  // transform coordinates to CommonTangentPlane
188  toCommonTangentPlane->transformPosAndErrors(*fs, *fs);
189  fittedStarList.push_back(fs);
190  mstar->setFittedStar(fs);
191  }
192  unMatchedCount++;
193  }
194  LOGLS_INFO(_log, "Unmatched objects: " << unMatchedCount);
195  } // end of loop on CcdImages
196 
197  // !!!!!!!!!!!!!!!!!
198  // TODO: DO WE REALLY NEED THIS???
199  // Why do we need to do this, instead of directly computing them in normalizeFittedStars?
200  // What makes the magnitudes special here?
201  // !!!!!!!!!!!!!!!!!
202  // assignMags();
203 }
std::unique_ptr< StarMatchList > listMatchCollect(const BaseStarList &list1, const BaseStarList &list2, const AstrometryTransform *guess, const double maxDist)
assembles star matches.
Definition: ListMatch.cc:569
#define LOGLS_DEBUG(logger, message)
Log a debug-level message using an iostream-based interface.
Definition: Log.h:608
T dynamic_pointer_cast(T... args)
T clear(T... args)
#define LOGLS_INFO(logger, message)
Log a info-level message using an iostream-based interface.
Definition: Log.h:628
T get(T... args)
BaseStarList & Fitted2Base(FittedStarList &This)
Definition: FittedStar.cc:64
BaseStarList & Measured2Base(MeasuredStarList &This)
Definition: MeasuredStar.cc:58
FittedStarList fittedStarList
Definition: Associations.h:58

◆ collectRefStars()

void lsst::jointcal::Associations::collectRefStars ( afw::table::SimpleCatalog refCat,
geom::Angle  matchCut,
std::string const &  fluxField,
float  refCoordinateErr,
bool  rejectBadFluxes = false 
)

Collect stars from an external reference catalog and associate them with fittedStars.

Parameters
refCatThe catalog of reference sources
[in]matchCutSeparation radius to match fitted and reference stars.
fluxFieldThe field name in refCat to get the flux from.
refCoordinateErrError on reference catalog coordinates [mas]. If not NaN, this overrides the coord_*_err values in the reference catalog itself. This value is divided by cos(dec) before being used for ra_err.
rejectBadFluxesReject reference sources with flux=NaN or 0 and/or fluxErr=NaN or 0. Typically false for astrometry and true for photometry.

Definition at line 205 of file Associations.cc.

207  {
208  if (refCat.size() == 0) {
209  throw(LSST_EXCEPT(pex::exceptions::InvalidParameterError,
210  " reference catalog is empty : stop here "));
211  }
212 
213  afw::table::CoordKey coordKey = refCat.getSchema()["coord"];
214  // Handle reference catalogs that don't have position errors.
215  afw::table::Key<float> raErrKey;
216  afw::table::Key<float> decErrKey;
217  if (std::isnan(refCoordinateErr)) {
218  raErrKey = refCat.getSchema()["coord_raErr"];
219  decErrKey = refCat.getSchema()["coord_decErr"];
220  }
221 
222  auto fluxKey = refCat.getSchema().find<double>(fluxField).key;
223  // Handle reference catalogs that don't have flux errors.
224  afw::table::Key<double> fluxErrKey;
225  try {
226  fluxErrKey = refCat.getSchema().find<double>(fluxField + "Err").key;
227  } catch (pex::exceptions::NotFoundError &) {
228  LOGLS_WARN(_log, "Flux error field ("
229  << fluxField << "Err"
230  << ") not found in reference catalog. Not using ref flux errors.");
231  }
232 
233  refStarList.clear();
234  for (size_t i = 0; i < refCat.size(); i++) {
235  auto const &record = refCat.get(i);
236 
237  auto coord = record->get(coordKey);
238  double flux = record->get(fluxKey);
239  double fluxErr;
240  if (fluxErrKey.isValid()) {
241  fluxErr = record->get(fluxErrKey);
242  } else {
244  }
245  double ra = lsst::geom::radToDeg(coord.getLongitude());
246  double dec = lsst::geom::radToDeg(coord.getLatitude());
247  auto star = std::make_shared<RefStar>(ra, dec, flux, fluxErr);
248 
249  if (std::isnan(refCoordinateErr)) {
250  star->vx = record->get(raErrKey);
251  star->vy = record->get(decErrKey);
252  } else {
253  // Compute and use the fake errors
254  star->vx = std::pow(refCoordinateErr / 1000. / 3600. / std::cos(coord.getLatitude()), 2);
255  star->vy = std::pow(refCoordinateErr / 1000. / 3600., 2);
256  }
257  // TODO: cook up a covariance as none of our current refcats have it
258  star->vxy = 0.;
259 
260  // Reject sources with non-finite fluxes and flux errors, and fluxErr=0 (which gives chi2=inf).
261  if (rejectBadFluxes && (!std::isfinite(flux) || !std::isfinite(fluxErr) || fluxErr <= 0)) continue;
262  refStarList.push_back(star);
263  }
264 
265  // project on CTP (i.e. RaDec2CTP), in degrees
266  AstrometryTransformLinear identity;
267  TanRaDecToPixel raDecToCommonTangentPlane(identity, _commonTangentPoint);
268 
269  associateRefStars(matchCut.asArcseconds(), &raDecToCommonTangentPlane);
270 }
#define LOGLS_WARN(logger, message)
Log a warn-level message using an iostream-based interface.
Definition: Log.h:648
double dec
Definition: Match.cc:41
constexpr double radToDeg(double x) noexcept
Definition: Angle.h:52
T push_back(T... args)
Key< U > key
Definition: Schema.cc:281
T isfinite(T... args)
T cos(T... args)
T clear(T... args)
constexpr double asArcseconds() const noexcept
Return an Angle&#39;s value in arcseconds.
Definition: Angle.h:178
T get(T... args)
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition: Exception.h:48
T pow(T... args)
T isnan(T... args)
T quiet_NaN(T... args)

◆ computeBoundingCircle()

lsst::sphgeom::Circle lsst::jointcal::Associations::computeBoundingCircle ( ) const

Return the bounding circle in on-sky (RA, Dec) coordinates containing all CcdImages.

Requires that computeCommonTangentPoint() be called first, so that sensor bounding boxes can be transformed into the common tangent plane.

Definition at line 92 of file Associations.cc.

92  {
93  // Compute the frame on the common tangent plane that contains all input images.
94  Frame tangentPlaneFrame;
95 
96  for (auto const &ccdImage : ccdImageList) {
97  Frame CTPFrame = ccdImage->getPixelToCommonTangentPlane()->apply(ccdImage->getImageFrame(), false);
98  if (tangentPlaneFrame.getArea() == 0)
99  tangentPlaneFrame = CTPFrame;
100  else
101  tangentPlaneFrame += CTPFrame;
102  }
103 
104  // Convert tangent plane coordinates to RaDec.
105  AstrometryTransformLinear identity;
106  TanPixelToRaDec commonTangentPlaneToRaDec(identity, _commonTangentPoint);
107  Frame raDecFrame = commonTangentPlaneToRaDec.apply(tangentPlaneFrame, false);
108 
110  points.reserve(4);
111  // raDecFrame is in on-sky (RA,Dec) degrees stored as an x/y box:
112  // the on-sky bounding box it represents is given by the corners of that box.
113  points.emplace_back(sphgeom::LonLat::fromDegrees(raDecFrame.xMin, raDecFrame.yMin));
114  points.emplace_back(sphgeom::LonLat::fromDegrees(raDecFrame.xMax, raDecFrame.yMin));
115  points.emplace_back(sphgeom::LonLat::fromDegrees(raDecFrame.xMin, raDecFrame.yMax));
116  points.emplace_back(sphgeom::LonLat::fromDegrees(raDecFrame.xMax, raDecFrame.yMax));
117 
119 }
Circle getBoundingCircle() const override
getBoundingCircle returns a bounding-circle for this region.
static LonLat fromDegrees(double lon, double lat)
Definition: LonLat.h:50
static ConvexPolygon convexHull(std::vector< UnitVector3d > const &points)
convexHull returns the convex hull of the given set of points if it exists and throws an exception ot...
Definition: ConvexPolygon.h:65
STL class.
T reserve(T... args)
T emplace_back(T... args)

◆ computeCommonTangentPoint()

void lsst::jointcal::Associations::computeCommonTangentPoint ( )

Sets a shared tangent point for all ccdImages, using the mean of the centers of all ccdImages.

Definition at line 76 of file Associations.cc.

76  {
78  centers.reserve(ccdImageList.size());
79  for (auto const &ccdImage : ccdImageList) {
80  centers.push_back(ccdImage->getBoresightRaDec());
81  }
82  auto commonTangentPoint = geom::averageSpherePoint(centers);
83  LOGLS_DEBUG(_log, "Using common tangent point: " << commonTangentPoint.getPosition(geom::degrees));
84  setCommonTangentPoint(commonTangentPoint.getPosition(geom::degrees));
85 }
T push_back(T... args)
#define LOGLS_DEBUG(logger, message)
Log a debug-level message using an iostream-based interface.
Definition: Log.h:608
AngleUnit constexpr degrees
constant with units of degrees
Definition: Angle.h:109
T size(T... args)
void setCommonTangentPoint(lsst::geom::Point2D const &commonTangentPoint)
Sets a shared tangent point for all ccdImages.
Definition: Associations.cc:87
STL class.
SpherePoint averageSpherePoint(std::vector< SpherePoint > const &coords)
Return the average of a list of coordinates.
Definition: SpherePoint.cc:235
T reserve(T... args)

◆ createCcdImage()

void lsst::jointcal::Associations::createCcdImage ( afw::table::SourceCatalog catalog,
std::shared_ptr< lsst::afw::geom::SkyWcs wcs,
std::shared_ptr< lsst::afw::image::VisitInfo visitInfo,
lsst::geom::Box2I const &  bbox,
std::string const &  filter,
std::shared_ptr< afw::image::PhotoCalib photoCalib,
std::shared_ptr< afw::cameraGeom::Detector detector,
int  visit,
int  ccd,
lsst::jointcal::JointcalControl const &  control 
)

Create a ccdImage from an exposure catalog and metadata, and add it to the list.

Parameters
[in]catalogThe extracted source catalog, selected for good astrometric sources.
[in]wcsThe exposure's original wcs
[in]visitInfoThe exposure's visitInfo object
[in]bboxThe bounding box of the exposure
[in]filterThe exposure's filter
[in]photoCalibThe exposure's photometric calibration
[in]detectorThe exposure's detector
[in]visitThe visit identifier
[in]ccdThe ccd identifier
[in]controlThe JointcalControl object

Definition at line 62 of file Associations.cc.

68  {
69  auto ccdImage = std::make_shared<CcdImage>(catalog, wcs, visitInfo, bbox, filter, photoCalib, detector,
70  visit, ccd, control.sourceFluxField);
72  LOGLS_DEBUG(_log, "Catalog " << ccdImage->getName() << " has " << ccdImage->getWholeCatalog().size()
73  << " objects.");
74 }
AmpInfoBoxKey bbox
Definition: Amplifier.cc:117
Key< int > visitInfo
Definition: Exposure.cc:70
Key< int > photoCalib
Definition: Exposure.cc:67
table::Key< table::Array< std::uint8_t > > wcs
Definition: SkyWcs.cc:71
T push_back(T... args)
#define LOGLS_DEBUG(logger, message)
Log a debug-level message using an iostream-based interface.
Definition: Log.h:608
table::Key< int > detector

◆ deprojectFittedStars()

void lsst::jointcal::Associations::deprojectFittedStars ( )

Sends back the fitted stars coordinates on the sky FittedStarsList::inTangentPlaneCoordinates keeps track of that.

Definition at line 394 of file Associations.cc.

394  {
395  // By default, Associations::fittedStarList is expressed on the Associations::commonTangentPlane.
396  // For AstrometryFit, we need it on the sky.
398  LOGLS_WARN(_log,
399  "DeprojectFittedStars: Fitted stars are already in sidereal coordinates, nothing done ");
400  return;
401  }
402 
403  TanPixelToRaDec ctp2Sky(AstrometryTransformLinear(), getCommonTangentPoint());
406 }
#define LOGLS_WARN(logger, message)
Log a warn-level message using an iostream-based interface.
Definition: Log.h:648
void applyTransform(const Operator &op)
enables to apply a geometrical transform if Star is Basestar or derives from it.
Definition: StarList.h:98
Point getCommonTangentPoint() const
can be used to project sidereal coordinates related to the image set on a plane.
Definition: Associations.h:105
FittedStarList fittedStarList
Definition: Associations.h:58

◆ fittedStarListSize()

size_t lsst::jointcal::Associations::fittedStarListSize ( )
inline

Definition at line 63 of file Associations.h.

63 { return fittedStarList.size(); }
T size(T... args)
FittedStarList fittedStarList
Definition: Associations.h:58

◆ getCcdImageList()

CcdImageList const& lsst::jointcal::Associations::getCcdImageList ( ) const
inline

Definition at line 171 of file Associations.h.

171 { return ccdImageList; }

◆ getCommonTangentPoint()

Point lsst::jointcal::Associations::getCommonTangentPoint ( ) const
inline

can be used to project sidereal coordinates related to the image set on a plane.

Definition at line 105 of file Associations.h.

105 { return _commonTangentPoint; }

◆ getNFilters()

unsigned lsst::jointcal::Associations::getNFilters ( ) const
inline

Number of different bands in the input image list. Not implemented so far.

Definition at line 174 of file Associations.h.

174 { return 1; }

◆ nCcdImagesValidForFit()

int lsst::jointcal::Associations::nCcdImagesValidForFit ( ) const

return the number of CcdImages with non-empty catalogs to-be-fit.

Definition at line 408 of file Associations.cc.

408  {
410  return item->getCatalogForFit().size() > 0;
411  });
412 }
T end(T... args)
T count_if(T... args)
T begin(T... args)

◆ nFittedStarsWithAssociatedRefStar()

size_t lsst::jointcal::Associations::nFittedStarsWithAssociatedRefStar ( ) const

Return the number of fittedStars that have an associated refStar.

Definition at line 414 of file Associations.cc.

414  {
415  size_t count = 0;
416  for (auto const &fittedStar : fittedStarList) {
417  if ((fittedStar != nullptr) & (fittedStar->getRefStar() != nullptr)) count++;
418  }
419  return count;
420 }
T count(T... args)
FittedStarList fittedStarList
Definition: Associations.h:58

◆ operator=() [1/2]

Associations& lsst::jointcal::Associations::operator= ( Associations const &  )
delete

◆ operator=() [2/2]

Associations& lsst::jointcal::Associations::operator= ( Associations &&  )
delete

◆ prepareFittedStars()

void lsst::jointcal::Associations::prepareFittedStars ( int  minMeasurements)

Set the color field of FittedStar 's from a colored catalog.

Prepare the fittedStar list by making quality cuts and normalizing measurements.

Parameters
[in]minMeasurementsThe minimum number of measuredStars for a FittedStar to be included.

Definition at line 298 of file Associations.cc.

298  {
299  selectFittedStars(minMeasurements);
300  normalizeFittedStars();
301 }

◆ refStarListSize()

size_t lsst::jointcal::Associations::refStarListSize ( )
inline

Definition at line 62 of file Associations.h.

62 { return refStarList.size(); }
T size(T... args)

◆ setCommonTangentPoint()

void lsst::jointcal::Associations::setCommonTangentPoint ( lsst::geom::Point2D const &  commonTangentPoint)

Sets a shared tangent point for all ccdImages.

Parameters
commonTangentPointThe common tangent point of all input images (decimal degrees).

Definition at line 87 of file Associations.cc.

87  {
88  _commonTangentPoint = Point(commonTangentPoint.getX(), commonTangentPoint.getY()); // a jointcal::Point
89  for (auto &ccdImage : ccdImageList) ccdImage->setCommonTangentPoint(_commonTangentPoint);
90 }

Member Data Documentation

◆ ccdImageList

CcdImageList lsst::jointcal::Associations::ccdImageList

Definition at line 56 of file Associations.h.

◆ fittedStarList

FittedStarList lsst::jointcal::Associations::fittedStarList

Definition at line 58 of file Associations.h.

◆ refStarList

RefStarList lsst::jointcal::Associations::refStarList

Definition at line 57 of file Associations.h.


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