LSST Applications g0f08755f38+9c285cab97,g1635faa6d4+13f3999e92,g1653933729+a8ce1bb630,g1a0ca8cf93+bf6eb00ceb,g28da252d5a+0829b12dee,g29321ee8c0+5700dc9eac,g2bbee38e9b+9634bc57db,g2bc492864f+9634bc57db,g2cdde0e794+c2c89b37c4,g3156d2b45e+41e33cbcdc,g347aa1857d+9634bc57db,g35bb328faa+a8ce1bb630,g3a166c0a6a+9634bc57db,g3e281a1b8c+9f2c4e2fc3,g414038480c+077ccc18e7,g41af890bb2+fde0dd39b6,g5fbc88fb19+17cd334064,g781aacb6e4+a8ce1bb630,g80478fca09+55a9465950,g82479be7b0+d730eedb7d,g858d7b2824+9c285cab97,g9125e01d80+a8ce1bb630,g9726552aa6+10f999ec6a,ga5288a1d22+2a84bb7594,gacf8899fa4+c69c5206e8,gae0086650b+a8ce1bb630,gb58c049af0+d64f4d3760,gc28159a63d+9634bc57db,gcf0d15dbbd+4b7d09cae4,gda3e153d99+9c285cab97,gda6a2b7d83+4b7d09cae4,gdaeeff99f8+1711a396fd,ge2409df99d+5e831397f4,ge79ae78c31+9634bc57db,gf0baf85859+147a0692ba,gf3967379c6+41c94011de,gf3fb38a9a8+8f07a9901b,gfb92a5be7c+9c285cab97,w.2024.46
LSST Data Management Base Package
Loading...
Searching...
No Matches
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.
 
 Associations (CcdImageList imageList, double epoch=0)
 Create an Associations object from a pre-built list of ccdImages.
 
 Associations (Associations const &)=delete
 No moves or copies: jointcal only ever needs one Associations object.
 
 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.
 
size_t getMaxMeasuredStars () const
 The number of MeasuredStars at the start of fitting, before any outliers are removed.
 
void createCcdImage (afw::table::SourceCatalog &catalog, const std::shared_ptr< lsst::afw::geom::SkyWcs > &wcs, const std::shared_ptr< lsst::afw::image::VisitInfo > &visitInfo, const lsst::geom::Box2I &bbox, const std::string &filter, const std::shared_ptr< afw::image::PhotoCalib > &photoCalib, const std::shared_ptr< afw::cameraGeom::Detector > &detector, int visit, int ccd, const lsst::jointcal::JointcalControl &control)
 Create a ccdImage from an exposure catalog and metadata, and add it to the list.
 
void addCcdImage (std::shared_ptr< CcdImage > const &ccdImage)
 Add a pre-constructed ccdImage to the ccdImageList.
 
void associateCatalogs (double matchCutInArcsec=0, bool useFittedList=false, bool enlargeFittedList=true)
 incrementaly builds a merged catalog of all image catalogs
 
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.
 
void deprojectFittedStars ()
 Sends back the fitted stars coordinates on the sky FittedStarsList::inTangentPlaneCoordinates keeps track of that.
 
void prepareFittedStars (int minMeasurements)
 Prepare the fittedStar list by making quality cuts and normalizing measurements.
 
void cleanFittedStars ()
 Remove FittedStars that have no measured stars; this can happen after outlier rejection.
 
CcdImageList const & getCcdImageList () const
 
unsigned getNFilters () const
 Number of different bands in the input image list. Not implemented so far.
 
lsst::sphgeom::Circle computeBoundingCircle () const
 Return the bounding circle in on-sky (RA, Dec) coordinates containing all CcdImages.
 
int nCcdImagesValidForFit () const
 return the number of CcdImages with non-empty catalogs to-be-fit.
 
size_t nFittedStarsWithAssociatedRefStar () const
 Return the number of fittedStars that have an associated refStar.
 
void setCommonTangentPoint (lsst::geom::Point2D const &commonTangentPoint)
 Shared tangent point for all ccdImages (decimal degrees).
 
Point getCommonTangentPoint () const
 Shared tangent point for all ccdImages (decimal degrees).
 
double getEpoch () const
 Common epoch of all of the ccdImages as a Julian Epoch Year (e.g.
 
void setEpoch (double epoch)
 Common epoch of all of the ccdImages as a Julian Epoch Year (e.g.
 

Public Attributes

CcdImageList ccdImageList
 
RefStarList refStarList
 
FittedStarList fittedStarList
 

Detailed Description

The class that implements the relations between MeasuredStar and FittedStar.

Definition at line 55 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 70 of file Associations.h.

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

◆ Associations() [2/4]

lsst::jointcal::Associations::Associations ( CcdImageList 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 84 of file Associations.h.

85 : ccdImageList(std::move(imageList)),
86 _commonTangentPoint(Point(std::numeric_limits<double>::quiet_NaN(),
88 _maxMeasuredStars(0),
89 _epoch(epoch) {}
T move(T... args)

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

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

◆ associateCatalogs()

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

incrementaly builds a merged catalog of all image catalogs

Definition at line 120 of file Associations.cc.

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

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

424 {
425 auto iter = fittedStarList.begin();
426 auto end = fittedStarList.end();
427 size_t count = 0;
428 while (iter != end) {
429 auto fittedStar = *iter;
430 if (fittedStar->getMeasurementCount() == 0) {
431 LOGLS_TRACE(_log, "Deleting FittedStar (has no measuredStars): " << *fittedStar);
432 iter = fittedStarList.erase(iter);
433 count++;
434 } else {
435 ++iter;
436 }
437 }
438 if (count > 0) {
439 LOGLS_INFO(_log, "Removed " << count << " fittedStars that had no associated measuredStar.");
440 }
441}
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 204 of file Associations.cc.

206 {
207 if (refCat.size() == 0) {
208 throw(LSST_EXCEPT(pex::exceptions::InvalidParameterError,
209 " reference catalog is empty : stop here "));
210 }
211
212 afw::table::CoordKey coordKey = refCat.getSchema()["coord"];
213 // Handle reference catalogs that don't have position errors.
214 afw::table::Key<float> raErrKey;
215 afw::table::Key<float> decErrKey;
216 if (std::isnan(refCoordinateErr)) {
217 raErrKey = refCat.getSchema()["coord_raErr"];
218 decErrKey = refCat.getSchema()["coord_decErr"];
219 }
220
221 auto fluxKey = refCat.getSchema().find<double>(fluxField).key;
222 // Handle reference catalogs that don't have flux errors.
223 afw::table::Key<double> fluxErrKey;
224 try {
225 fluxErrKey = refCat.getSchema().find<double>(fluxField + "Err").key;
226 } catch (pex::exceptions::NotFoundError &) {
227 LOGLS_WARN(_log, "Flux error field ("
228 << fluxField << "Err"
229 << ") not found in reference catalog. Not using ref flux errors.");
230 }
231
232 // Handle reference catalogs that don't have proper motion & error
233 afw::table::Key<geom::Angle> pmRaKey, pmDecKey;
234 afw::table::Key<float> pmRaErrKey, pmDecErrKey, pmRaDecCovKey;
235 try {
236 pmRaKey = refCat.getSchema().find<geom::Angle>("pm_ra").key;
237 pmDecKey = refCat.getSchema().find<geom::Angle>("pm_dec").key;
238 } catch (pex::exceptions::NotFoundError &ex) {
239 LOGLS_WARN(_log, "Not loading proper motions: (pm_ra,pm_dec) fields not found in reference catalog.");
240 } catch (pex::exceptions::TypeError &ex) {
241 LOGLS_WARN(_log, "Not loading proper motions: RA/Dec proper motion values must be `geom:Angle`: "
242 << ex.what());
243 }
244 try {
245 pmRaErrKey = refCat.getSchema().find<float>("pm_raErr").key;
246 pmDecErrKey = refCat.getSchema().find<float>("pm_decErr").key;
247 } catch (pex::exceptions::NotFoundError &ex) {
248 LOGLS_WARN(_log, "Not loading proper motions: error fields not available: " << ex.what());
249 }
250
251 // TODO: we aren't getting covariances from Gaia yet, so maybe ignore this for now?
252 try {
253 pmRaDecCovKey = refCat.getSchema().find<float>("pm_ra_pm_dec_Cov").key;
254 } catch (pex::exceptions::NotFoundError &ex) {
255 LOGLS_WARN(_log, "No ra/dec proper motion covariances in refcat: " << ex.what());
256 }
257
259 for (size_t i = 0; i < refCat.size(); i++) {
260 auto const &record = refCat.get(i);
261
262 auto coord = record->get(coordKey);
263 double flux = record->get(fluxKey);
264 double fluxErr;
265 if (fluxErrKey.isValid()) {
266 fluxErr = record->get(fluxErrKey);
267 } else {
269 }
270 double ra = lsst::geom::radToDeg(coord.getLongitude());
271 double dec = lsst::geom::radToDeg(coord.getLatitude());
272 auto star = std::make_shared<RefStar>(ra, dec, flux, fluxErr);
273
274 if (std::isnan(refCoordinateErr)) {
275 // refcat errors are unitless but stored as radians: convert to deg**2
276 star->vx = std::pow(lsst::geom::radToDeg(record->get(raErrKey)), 2);
277 star->vy = std::pow(lsst::geom::radToDeg(record->get(decErrKey)), 2);
278 } else {
279 // Convert the fake errors from mas to deg**2
280 star->vx = std::pow(refCoordinateErr / 1000. / 3600. / std::cos(coord.getLatitude()), 2);
281 star->vy = std::pow(refCoordinateErr / 1000. / 3600., 2);
282 }
283
284 if (pmRaKey.isValid()) {
285 if (pmRaDecCovKey.isValid()) {
286 star->setProperMotion(std::make_unique<ProperMotion const>(
287 record->get(pmRaKey).asRadians(), record->get(pmDecKey).asRadians(),
288 record->get(pmRaErrKey), record->get(pmDecErrKey), record->get(pmRaDecCovKey)));
289 } else {
290 star->setProperMotion(std::make_unique<ProperMotion const>(
291 record->get(pmRaKey).asRadians(), record->get(pmDecKey).asRadians(),
292 record->get(pmRaErrKey), record->get(pmDecErrKey)));
293 }
294 }
295
296 // TODO: cook up a covariance as none of our current refcats have it
297 star->vxy = 0.;
298
299 // Reject sources with non-finite fluxes and flux errors, and fluxErr=0 (which gives chi2=inf).
300 if (rejectBadFluxes && (!std::isfinite(flux) || !std::isfinite(fluxErr) || fluxErr <= 0)) continue;
302 }
303
304 // project on CTP (i.e. RaDec2CTP), in degrees
305 AstrometryTransformLinear identity;
306 TanRaDecToPixel raDecToCommonTangentPlane(identity, _commonTangentPoint);
307
308 associateRefStars(matchCut.asArcseconds(), &raDecToCommonTangentPlane);
309}
#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:128
constexpr double asArcseconds() const noexcept
Return an Angle's value in arcseconds.
Definition Angle.h:185
T cos(T... args)
T isfinite(T... args)
T isnan(T... args)
constexpr double radToDeg(double x) noexcept
Definition Angle.h:53
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 LOGLS_INFO(_log, "Computed tangent plane box for data (degrees): " << raDecFrame);
112 points.emplace_back(sphgeom::LonLat::fromDegrees(raDecFrame.xMin, raDecFrame.yMin));
113 points.emplace_back(sphgeom::LonLat::fromDegrees(raDecFrame.xMax, raDecFrame.yMin));
114 points.emplace_back(sphgeom::LonLat::fromDegrees(raDecFrame.xMin, raDecFrame.yMax));
115 points.emplace_back(sphgeom::LonLat::fromDegrees(raDecFrame.xMax, raDecFrame.yMax));
116
118}
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...
Circle getBoundingCircle() const override
getBoundingCircle returns a bounding-circle for this region.
static LonLat fromDegrees(double lon, double lat)
Definition LonLat.h:57
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).
SpherePoint averageSpherePoint(std::vector< SpherePoint > const &coords)
Return the average of a list of coordinates.
AngleUnit constexpr degrees
constant with units of degrees
Definition Angle.h:110
T size(T... args)

◆ createCcdImage()

void lsst::jointcal::Associations::createCcdImage ( afw::table::SourceCatalog & catalog,
const std::shared_ptr< lsst::afw::geom::SkyWcs > & wcs,
const std::shared_ptr< lsst::afw::image::VisitInfo > & visitInfo,
const lsst::geom::Box2I & bbox,
const std::string & filter,
const std::shared_ptr< afw::image::PhotoCalib > & photoCalib,
const std::shared_ptr< afw::cameraGeom::Detector > & detector,
int visit,
int ccd,
const lsst::jointcal::JointcalControl & 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
std::string sourceFluxField
"name of flux field in source catalog" ;
Key< int > visitInfo
Definition Exposure.cc:70

◆ deprojectFittedStars()

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

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

Definition at line 454 of file Associations.cc.

454 {
455 // By default, Associations::fittedStarList is expressed on the Associations::commonTangentPlane.
456 // For AstrometryFit, we need it on the sky.
458 LOGLS_WARN(_log,
459 "DeprojectFittedStars: Fitted stars are already in sidereal coordinates, nothing done ");
460 return;
461 }
462
463 TanPixelToRaDec ctp2Sky(AstrometryTransformLinear(), getCommonTangentPoint());
466}
Point getCommonTangentPoint() const
Shared tangent point for all ccdImages (decimal degrees).
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 64 of file Associations.h.

64{ return fittedStarList.size(); }

◆ getCcdImageList()

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

Definition at line 189 of file Associations.h.

189{ 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 109 of file Associations.h.

109{ 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 119 of file Associations.h.

119{ 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 113 of file Associations.h.

113{ 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 192 of file Associations.h.

192{ return 1; }

◆ nCcdImagesValidForFit()

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

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

Definition at line 468 of file Associations.cc.

468 {
470 return item->getCatalogForFit().size() > 0;
471 });
472}

◆ nFittedStarsWithAssociatedRefStar()

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

Return the number of fittedStars that have an associated refStar.

Definition at line 474 of file Associations.cc.

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

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

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

◆ refStarListSize()

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

Definition at line 63 of file Associations.h.

63{ 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 120 of file Associations.h.

120{ _epoch = epoch; }

Member Data Documentation

◆ ccdImageList

CcdImageList lsst::jointcal::Associations::ccdImageList

Definition at line 57 of file Associations.h.

◆ fittedStarList

FittedStarList lsst::jointcal::Associations::fittedStarList

Definition at line 59 of file Associations.h.

◆ refStarList

RefStarList lsst::jointcal::Associations::refStarList

Definition at line 58 of file Associations.h.


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