6 #include "Eigen/Sparse" 18 LOG_LOGGER _log =
LOG_GET(
"jointcal.PhotometryFit");
24 void PhotometryFit::leastSquareDerivativesMeasurement(CcdImage
const &
ccdImage, TripletList &tripletList,
25 Eigen::VectorXd &grad,
26 MeasuredStarList
const *measuredStarList)
const {
35 if (measuredStarList) assert(&(measuredStarList->front()->getCcdImage()) == &ccdImage);
37 unsigned nparModel = (_fittingModel) ? _photometryModel->getNpar(ccdImage) : 0;
38 unsigned nparFlux = (_fittingFluxes) ? 1 : 0;
39 unsigned nparTotal = nparModel + nparFlux;
41 if (_fittingModel) _photometryModel->getMappingIndices(ccdImage, indices);
43 Eigen::VectorXd H(nparTotal);
45 unsigned kTriplets = tripletList.getNextFreeIndex();
46 const MeasuredStarList &catalog = (measuredStarList) ? *measuredStarList : ccdImage.getCatalogForFit();
48 for (
auto const &measuredStar : catalog) {
49 if (!measuredStar->isValid())
continue;
52 TweakPhotomMeasurementErrors(inPos, *measuredStar, _fluxError);
56 double residual = _photometryModel->transform(ccdImage, *measuredStar, measuredStar->getInstFlux()) -
57 measuredStar->getFittedStar()->getFlux();
60 1.0 / _photometryModel->transform(ccdImage, *measuredStar, measuredStar->getInstFluxErr());
61 double W =
std::pow(inverseSigma, 2);
64 _photometryModel->computeParameterDerivatives(*measuredStar, ccdImage, H);
65 for (
unsigned k = 0; k < indices.size(); k++) {
66 unsigned l = indices[k];
67 tripletList.addTriplet(l, kTriplets, H[k] * inverseSigma);
68 grad[l] += H[k] * W * residual;
72 unsigned index = measuredStar->getFittedStar()->getIndexInMatrix();
74 tripletList.addTriplet(index, kTriplets, -1.0 * inverseSigma);
75 grad[index] += -1.0 * W * residual;
80 tripletList.setNextFreeIndex(kTriplets);
83 void PhotometryFit::leastSquareDerivativesReference(FittedStarList
const &fittedStarList,
84 TripletList &tripletList, Eigen::VectorXd &grad)
const {
91 if (!_fittingFluxes)
return;
95 unsigned kTriplets = tripletList.getNextFreeIndex();
97 for (
auto const &fittedStar : fittedStarList) {
98 auto refStar = fittedStar->getRefStar();
99 if (refStar ==
nullptr)
continue;
106 double inverseSigma = 1.0 / refStar->getFluxErr();
108 double residual = fittedStar->getFlux() - refStar->getFlux();
110 unsigned index = fittedStar->getIndexInMatrix();
112 tripletList.addTriplet(index, kTriplets, 1.0 * inverseSigma);
113 grad(index) += 1.0 *
std::pow(inverseSigma, 2) * residual;
116 tripletList.setNextFreeIndex(kTriplets);
119 void PhotometryFit::accumulateStatImageList(
CcdImageList const &ccdImageList, Chi2Accumulator &accum)
const {
124 for (
auto const &ccdImage : ccdImageList) {
125 auto &catalog = ccdImage->getCatalogForFit();
127 for (
auto const &measuredStar : catalog) {
128 if (!measuredStar->isValid())
continue;
130 _photometryModel->transform(*ccdImage, *measuredStar, measuredStar->getInstFluxErr());
132 TweakPhotomMeasurementErrors(inPos, measuredStar, _fluxError);
135 _photometryModel->transform(*ccdImage, *measuredStar, measuredStar->getInstFlux()) -
136 measuredStar->getFittedStar()->getFlux();
138 double chi2Val =
std::pow(residual / sigma, 2);
139 accum.addEntry(chi2Val, 1, measuredStar);
144 void PhotometryFit::accumulateStatRefStars(Chi2Accumulator &accum)
const {
150 FittedStarList &fittedStarList =
_associations->fittedStarList;
151 for (
auto const &fittedStar : fittedStarList) {
152 auto refStar = fittedStar->getRefStar();
153 if (refStar ==
nullptr)
continue;
154 double chi2 =
std::pow(((fittedStar->getFlux() - refStar->getFlux()) / refStar->getFluxErr()), 2);
155 accum.addEntry(chi2, 1, fittedStar);
162 void PhotometryFit::getIndicesOfMeasuredStar(MeasuredStar
const &measuredStar,
166 _photometryModel->getMappingIndices(measuredStar.getCcdImage(), indices);
168 if (_fittingFluxes) {
169 auto fs = measuredStar.getFittedStar();
170 unsigned fsIndex = fs->getIndexInMatrix();
177 LOGLS_INFO(_log,
"assignIndices: now fitting: " << whatToFit);
179 _fittingFluxes = (
_whatToFit.
find(
"Fluxes") != std::string::npos);
182 _nParModel = (_fittingModel) ? _photometryModel->assignIndices(whatToFit, 0) : 0;
183 unsigned ipar = _nParModel;
185 if (_fittingFluxes) {
191 fittedStar->setIndexInMatrix(ipar);
195 _nParFluxes = ipar - _nParModel;
198 "nParameters total: " <<
_nParTot <<
" model: " << _nParModel <<
" fluxes: " << _nParFluxes);
204 "PhotometryFit::offsetParams : the provided vector length is not compatible with " 205 "the current whatToFit setting");
206 if (_fittingModel) _photometryModel->offsetParams(delta);
208 if (_fittingFluxes) {
213 unsigned index = fittedStar->getIndexInMatrix();
214 fittedStar->getFlux() -= delta(index);
225 ofile <<
"#xccd: coordinate in CCD" 231 <<
"instFlux : measured instrument flux" 233 <<
"instFluxError : measured instrument flux error" 235 <<
"flux : measured flux" 237 <<
"fluxError : measured flux error" 239 <<
"transformedFlux:" 241 <<
"transformedFluxErr:" 243 <<
"fflux : fitted flux" 245 <<
"mjd: Julian date of the measurement" 249 <<
"fsindex: some unique index of the object" 251 <<
"ra: pos of fitted star" 253 <<
"dec: pos of fitted star" 255 <<
"chi2: contribution to Chi2 (1 dof)" 257 <<
"nm: number of measurements of this FittedStar" 259 <<
"chip: chip number" 265 for (
auto const &i : ccdImageList) {
268 for (
auto const &is : cat) {
270 if (!measuredStar.
isValid())
continue;
273 tweakPhotomMeasurementErrors(inPos, measuredStar, _fluxError);
275 double flux = _photometryModel->transform(ccdImage, measuredStar, measuredStar.
getInstFlux());
277 _photometryModel->transform(ccdImage, measuredStar, measuredStar.
getInstFluxErr());
278 double jd = ccdImage.
getMjd();
280 double residual = flux - fs->getFlux();
281 double chi2Val =
std::pow(residual / sigma, 2);
282 ofile << measuredStar.
x <<
"\t" << measuredStar.
y <<
"\t" << fs->getMag() <<
"\t" 284 << measuredStar.
getFlux() <<
"\t" << measuredStar.
getFluxErr() <<
"\t" << flux <<
"\t" 285 << fluxErr <<
"\t" << fs->getFlux() <<
"\t" << jd <<
"\t" << fs->color <<
"\t" 286 << fs->getIndexInMatrix() <<
"\t" << fs->x <<
"\t" << fs->y <<
"\t" << chi2Val <<
"\t" 287 << fs->getMeasurementCount() <<
"\t" << ccdImage.
getCcdId() <<
"\t" << ccdImage.
getVisit()
VisitIdType getVisit() const
returns visit ID
double getMjd() const
Julian Date.
CcdIdType getCcdId() const
returns ccd ID
double getFluxErr() const
A list of MeasuredStar. They are usually filled in Associations::AddImage.
afw::table::Key< double > sigma
#define LOGLS_INFO(logger, message)
void offsetParams(Eigen::VectorXd const &delta) override
Offset the parameters by the requested quantities.
Class for a simple mapping implementing a generic Gtransfo.
void saveResultTuples(std::string const &tupleName) const override
Save the full chi2 term per star that was used in the minimization, for debugging.
double getInstFlux() const
double getInstFluxErr() const
objects measured on actual images.
std::shared_ptr< Associations > _associations
void assignIndices(std::string const &whatToFit) override
Set parameters to fit and assign indices in the big matrix.
#define LSST_EXCEPT(type,...)
std::shared_ptr< const FittedStar > getFittedStar() const
MeasuredStarList const & getCatalogForFit() const
Gets the catalog to be used for fitting, which may have been cleaned-up.
#define LOGLS_DEBUG(logger, message)
std::list< std::shared_ptr< CcdImage > > CcdImageList
Handler of an actual image from a single CCD.
bool isValid() const
Fits may use that to discard outliers.