LSST Applications  21.0.0-172-gfb10e10a+18fedfabac,22.0.0+297cba6710,22.0.0+80564b0ff1,22.0.0+8d77f4f51a,22.0.0+a28f4c53b1,22.0.0+dcf3732eb2,22.0.1-1-g7d6de66+2a20fdde0d,22.0.1-1-g8e32f31+297cba6710,22.0.1-1-geca5380+7fa3b7d9b6,22.0.1-12-g44dc1dc+2a20fdde0d,22.0.1-15-g6a90155+515f58c32b,22.0.1-16-g9282f48+790f5f2caa,22.0.1-2-g92698f7+dcf3732eb2,22.0.1-2-ga9b0f51+7fa3b7d9b6,22.0.1-2-gd1925c9+bf4f0e694f,22.0.1-24-g1ad7a390+a9625a72a8,22.0.1-25-g5bf6245+3ad8ecd50b,22.0.1-25-gb120d7b+8b5510f75f,22.0.1-27-g97737f7+2a20fdde0d,22.0.1-32-gf62ce7b1+aa4237961e,22.0.1-4-g0b3f228+2a20fdde0d,22.0.1-4-g243d05b+871c1b8305,22.0.1-4-g3a563be+32dcf1063f,22.0.1-4-g44f2e3d+9e4ab0f4fa,22.0.1-42-gca6935d93+ba5e5ca3eb,22.0.1-5-g15c806e+85460ae5f3,22.0.1-5-g58711c4+611d128589,22.0.1-5-g75bb458+99c117b92f,22.0.1-6-g1c63a23+7fa3b7d9b6,22.0.1-6-g50866e6+84ff5a128b,22.0.1-6-g8d3140d+720564cf76,22.0.1-6-gd805d02+cc5644f571,22.0.1-8-ge5750ce+85460ae5f3,master-g6e05de7fdc+babf819c66,master-g99da0e417a+8d77f4f51a,w.2021.48
LSST Data Management Base Package
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, double epoch=0)
 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...
 
size_t getMaxMeasuredStars () const
 The number of MeasuredStars at the start of fitting, before any outliers are removed. 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)
 Prepare the fittedStar list by making quality cuts and normalizing measurements. More...
 
void cleanFittedStars ()
 Remove FittedStars that have no measured stars; this can happen after outlier rejection. 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...
 
void setCommonTangentPoint (lsst::geom::Point2D const &commonTangentPoint)
 Shared tangent point for all ccdImages (decimal degrees). More...
 
Point getCommonTangentPoint () const
 Shared tangent point for all ccdImages (decimal degrees). More...
 
double getEpoch () const
 Common epoch of all of the ccdImages as a Julian Epoch Year (e.g. More...
 
void setEpoch (double epoch)
 Common epoch of all of the ccdImages as a Julian Epoch Year (e.g. 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.

70  : _commonTangentPoint(Point(std::numeric_limits<double>::quiet_NaN(),
72  _maxMeasuredStars(0) {}

◆ Associations() [2/4]

lsst::jointcal::Associations::Associations ( CcdImageList const &  imageList,
double  epoch = 0 
)
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.
epochThe julian epoch year to which all proper motion corrections should be made.

Definition at line 83 of file Associations.h.

84  : ccdImageList(imageList),
85  _commonTangentPoint(Point(std::numeric_limits<double>::quiet_NaN(),
87  _maxMeasuredStars(0),
88  _epoch(epoch) {}

◆ 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 145 of file Associations.h.

145 { 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 119 of file Associations.cc.

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

◆ cleanFittedStars()

void lsst::jointcal::Associations::cleanFittedStars ( )

Remove FittedStars that have no measured stars; this can happen after outlier rejection.

Use this to perform e.g. position minimization with outlier rejection after model minimization has removed measuredStar outliers, to prevent the matrix from becoming non-positive definite. Once this has been called, prepareFittedStars() has to be called again if the full measuredStar->FittedStar relationship needs to be reconstructed.

Definition at line 423 of file Associations.cc.

423  {
424  auto iter = fittedStarList.begin();
425  auto end = fittedStarList.end();
426  size_t count = 0;
427  while (iter != end) {
428  auto fittedStar = *iter;
429  if (fittedStar->getMeasurementCount() == 0) {
430  LOGLS_TRACE(_log, "Deleting FittedStar (has no measuredStars): " << *fittedStar);
432  count++;
433  } else {
434  ++iter;
435  }
436  }
437  if (count > 0) {
438  LOGLS_INFO(_log, "Removed " << count << " fittedStars that had no associated measuredStar.");
439  }
440 }
int end
#define LOGLS_INFO(logger, message)
Log a info-level message using an iostream-based interface.
Definition: Log.h:639
#define LOGLS_TRACE(logger, message)
Log a trace-level message using an iostream-based interface.
Definition: Log.h:599
T begin(T... args)
T count(T... args)
T end(T... args)
T erase(T... args)

◆ 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 203 of file Associations.cc.

205  {
206  if (refCat.size() == 0) {
207  throw(LSST_EXCEPT(pex::exceptions::InvalidParameterError,
208  " reference catalog is empty : stop here "));
209  }
210 
211  afw::table::CoordKey coordKey = refCat.getSchema()["coord"];
212  // Handle reference catalogs that don't have position errors.
213  afw::table::Key<float> raErrKey;
214  afw::table::Key<float> decErrKey;
215  if (std::isnan(refCoordinateErr)) {
216  raErrKey = refCat.getSchema()["coord_raErr"];
217  decErrKey = refCat.getSchema()["coord_decErr"];
218  }
219 
220  auto fluxKey = refCat.getSchema().find<double>(fluxField).key;
221  // Handle reference catalogs that don't have flux errors.
222  afw::table::Key<double> fluxErrKey;
223  try {
224  fluxErrKey = refCat.getSchema().find<double>(fluxField + "Err").key;
225  } catch (pex::exceptions::NotFoundError &) {
226  LOGLS_WARN(_log, "Flux error field ("
227  << fluxField << "Err"
228  << ") not found in reference catalog. Not using ref flux errors.");
229  }
230 
231  // Handle reference catalogs that don't have proper motion & error
232  afw::table::Key<geom::Angle> pmRaKey, pmDecKey;
233  afw::table::Key<float> pmRaErrKey, pmDecErrKey, pmRaDecCovKey;
234  try {
235  pmRaKey = refCat.getSchema().find<geom::Angle>("pm_ra").key;
236  pmDecKey = refCat.getSchema().find<geom::Angle>("pm_dec").key;
237  } catch (pex::exceptions::NotFoundError &ex) {
238  LOGLS_WARN(_log, "Not loading proper motions: (pm_ra,pm_dec) fields not found in reference catalog.");
239  } catch (pex::exceptions::TypeError &ex) {
240  LOGLS_WARN(_log, "Not loading proper motions: RA/Dec proper motion values must be `geom:Angle`: "
241  << ex.what());
242  }
243  try {
244  pmRaErrKey = refCat.getSchema().find<float>("pm_raErr").key;
245  pmDecErrKey = refCat.getSchema().find<float>("pm_decErr").key;
246  } catch (pex::exceptions::NotFoundError &ex) {
247  LOGLS_WARN(_log, "Not loading proper motions: error fields not available: " << ex.what());
248  }
249 
250  // TODO: we aren't getting covariances from Gaia yet, so maybe ignore this for now?
251  try {
252  pmRaDecCovKey = refCat.getSchema().find<float>("pm_ra_Dec_Cov").key;
253  } catch (pex::exceptions::NotFoundError &ex) {
254  LOGLS_WARN(_log, "No ra/dec proper motion covariances in refcat: " << ex.what());
255  }
256 
257  refStarList.clear();
258  for (size_t i = 0; i < refCat.size(); i++) {
259  auto const &record = refCat.get(i);
260 
261  auto coord = record->get(coordKey);
262  double flux = record->get(fluxKey);
263  double fluxErr;
264  if (fluxErrKey.isValid()) {
265  fluxErr = record->get(fluxErrKey);
266  } else {
268  }
269  double ra = lsst::geom::radToDeg(coord.getLongitude());
270  double dec = lsst::geom::radToDeg(coord.getLatitude());
271  auto star = std::make_shared<RefStar>(ra, dec, flux, fluxErr);
272 
273  if (std::isnan(refCoordinateErr)) {
274  // refcat errors are unitless but stored as radians: convert to deg**2
275  star->vx = std::pow(lsst::geom::radToDeg(record->get(raErrKey)), 2);
276  star->vy = std::pow(lsst::geom::radToDeg(record->get(decErrKey)), 2);
277  } else {
278  // Convert the fake errors from mas to deg**2
279  star->vx = std::pow(refCoordinateErr / 1000. / 3600. / std::cos(coord.getLatitude()), 2);
280  star->vy = std::pow(refCoordinateErr / 1000. / 3600., 2);
281  }
282 
283  if (pmRaKey.isValid()) {
284  if (pmRaDecCovKey.isValid()) {
285  star->setProperMotion(std::make_unique<ProperMotion const>(
286  record->get(pmRaKey).asRadians(), record->get(pmDecKey).asRadians(),
287  record->get(pmRaErrKey), record->get(pmDecErrKey), record->get(pmRaDecCovKey)));
288  } else {
289  star->setProperMotion(std::make_unique<ProperMotion const>(
290  record->get(pmRaKey).asRadians(), record->get(pmDecKey).asRadians(),
291  record->get(pmRaErrKey), record->get(pmDecErrKey)));
292  }
293  }
294 
295  // TODO: cook up a covariance as none of our current refcats have it
296  star->vxy = 0.;
297 
298  // Reject sources with non-finite fluxes and flux errors, and fluxErr=0 (which gives chi2=inf).
299  if (rejectBadFluxes && (!std::isfinite(flux) || !std::isfinite(fluxErr) || fluxErr <= 0)) continue;
300  refStarList.push_back(star);
301  }
302 
303  // project on CTP (i.e. RaDec2CTP), in degrees
304  AstrometryTransformLinear identity;
305  TanRaDecToPixel raDecToCommonTangentPlane(identity, _commonTangentPoint);
306 
307  associateRefStars(matchCut.asArcseconds(), &raDecToCommonTangentPlane);
308 }
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition: Exception.h:48
#define LOGLS_WARN(logger, message)
Log a warn-level message using an iostream-based interface.
Definition: Log.h:659
double dec
Definition: Match.cc:41
A class representing an angle.
Definition: Angle.h:127
constexpr double asArcseconds() const noexcept
Return an Angle's value in arcseconds.
Definition: Angle.h:178
T cos(T... args)
T isfinite(T... args)
T isnan(T... args)
constexpr double radToDeg(double x) noexcept
Definition: Angle.h:52
T pow(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 90 of file Associations.cc.

90  {
91  // Compute the frame on the common tangent plane that contains all input images.
92  Frame tangentPlaneFrame;
93 
94  for (auto const &ccdImage : ccdImageList) {
95  Frame CTPFrame = ccdImage->getPixelToCommonTangentPlane()->apply(ccdImage->getImageFrame(), false);
96  if (tangentPlaneFrame.getArea() == 0)
97  tangentPlaneFrame = CTPFrame;
98  else
99  tangentPlaneFrame += CTPFrame;
100  }
101 
102  // Convert tangent plane coordinates to RaDec.
103  AstrometryTransformLinear identity;
104  TanPixelToRaDec commonTangentPlaneToRaDec(identity, _commonTangentPoint);
105  Frame raDecFrame = commonTangentPlaneToRaDec.apply(tangentPlaneFrame, false);
106 
108  points.reserve(4);
109  // raDecFrame is in on-sky (RA,Dec) degrees stored as an x/y box:
110  // the on-sky bounding box it represents is given by the corners of that box.
111  points.emplace_back(sphgeom::LonLat::fromDegrees(raDecFrame.xMin, raDecFrame.yMin));
112  points.emplace_back(sphgeom::LonLat::fromDegrees(raDecFrame.xMax, raDecFrame.yMin));
113  points.emplace_back(sphgeom::LonLat::fromDegrees(raDecFrame.xMin, raDecFrame.yMax));
114  points.emplace_back(sphgeom::LonLat::fromDegrees(raDecFrame.xMax, raDecFrame.yMax));
115 
117 }
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
Circle getBoundingCircle() const override
getBoundingCircle returns a bounding-circle for this region.
static LonLat fromDegrees(double lon, double lat)
Definition: LonLat.h:50
T emplace_back(T... args)
T reserve(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 74 of file Associations.cc.

74  {
76  centers.reserve(ccdImageList.size());
77  for (auto const &ccdImage : ccdImageList) {
78  centers.push_back(ccdImage->getBoresightRaDec());
79  }
80  auto commonTangentPoint = geom::averageSpherePoint(centers);
81  LOGLS_DEBUG(_log, "Using common tangent point: " << commonTangentPoint.getPosition(geom::degrees));
82  setCommonTangentPoint(commonTangentPoint.getPosition(geom::degrees));
83 }
void setCommonTangentPoint(lsst::geom::Point2D const &commonTangentPoint)
Shared tangent point for all ccdImages (decimal degrees).
Definition: Associations.cc:85
SpherePoint averageSpherePoint(std::vector< SpherePoint > const &coords)
Return the average of a list of coordinates.
Definition: SpherePoint.cc:235
constexpr AngleUnit degrees
constant with units of degrees
Definition: Angle.h:109
T size(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);
71  ccdImageList.push_back(ccdImage);
72 }
AmpInfoBoxKey bbox
Definition: Amplifier.cc:117
table::Key< int > detector
table::Key< table::Array< std::uint8_t > > wcs
Definition: SkyWcs.cc:66
Key< int > visitInfo
Definition: Exposure.cc:70
Key< int > photoCalib
Definition: Exposure.cc:67

◆ deprojectFittedStars()

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

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

Definition at line 453 of file Associations.cc.

453  {
454  // By default, Associations::fittedStarList is expressed on the Associations::commonTangentPlane.
455  // For AstrometryFit, we need it on the sky.
457  LOGLS_WARN(_log,
458  "DeprojectFittedStars: Fitted stars are already in sidereal coordinates, nothing done ");
459  return;
460  }
461 
462  TanPixelToRaDec ctp2Sky(AstrometryTransformLinear(), getCommonTangentPoint());
465 }
Point getCommonTangentPoint() const
Shared tangent point for all ccdImages (decimal degrees).
Definition: Associations.h:108
void applyTransform(const Operator &op)
enables to apply a geometrical transform if Star is Basestar or derives from it.
Definition: StarList.h:98

◆ fittedStarListSize()

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

Definition at line 63 of file Associations.h.

63 { return fittedStarList.size(); }

◆ getCcdImageList()

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

Definition at line 188 of file Associations.h.

188 { return ccdImageList; }

◆ getCommonTangentPoint()

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

Shared tangent point for all ccdImages (decimal degrees).

Used to project sidereal coordinates related to the images onto a single plane.

Definition at line 108 of file Associations.h.

108 { return _commonTangentPoint; }

◆ getEpoch()

double lsst::jointcal::Associations::getEpoch ( ) const
inline

Common epoch of all of the ccdImages as a Julian Epoch Year (e.g.

2000.0 for J2000).

Definition at line 118 of file Associations.h.

118 { return _epoch; }

◆ getMaxMeasuredStars()

size_t lsst::jointcal::Associations::getMaxMeasuredStars ( ) const
inline

The number of MeasuredStars at the start of fitting, before any outliers are removed.

Definition at line 112 of file Associations.h.

112 { return _maxMeasuredStars; }

◆ getNFilters()

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

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

Definition at line 191 of file Associations.h.

191 { return 1; }

◆ nCcdImagesValidForFit()

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

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

Definition at line 467 of file Associations.cc.

467  {
469  return item->getCatalogForFit().size() > 0;
470  });
471 }

◆ nFittedStarsWithAssociatedRefStar()

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

Return the number of fittedStars that have an associated refStar.

Definition at line 473 of file Associations.cc.

473  {
474  size_t count = 0;
475  for (auto const &fittedStar : fittedStarList) {
476  if ((fittedStar != nullptr) & (fittedStar->getRefStar() != nullptr)) count++;
477  }
478  return count;
479 }

◆ operator=() [1/2]

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

◆ operator=() [2/2]

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

◆ prepareFittedStars()

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

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 336 of file Associations.cc.

336  {
337  selectFittedStars(minMeasurements);
338  normalizeFittedStars();
339 }

◆ refStarListSize()

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

Definition at line 62 of file Associations.h.

62 { return refStarList.size(); }

◆ setCommonTangentPoint()

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

Shared tangent point for all ccdImages (decimal degrees).

Used to project sidereal coordinates related to the images onto a single plane.

Definition at line 85 of file Associations.cc.

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

◆ setEpoch()

void lsst::jointcal::Associations::setEpoch ( double  epoch)
inline

Common epoch of all of the ccdImages as a Julian Epoch Year (e.g.

2000.0 for J2000).

Definition at line 119 of file Associations.h.

119 { _epoch = epoch; }

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: