31 #include "Eigen/Sparse" 
   45 void PhotometryFit::leastSquareDerivativesMeasurement(CcdImage 
const &ccdImage, TripletList &tripletList,
 
   46                                                       Eigen::VectorXd &grad,
 
   47                                                       MeasuredStarList 
const *measuredStarList)
 const {
 
   56     if (measuredStarList) assert(&(measuredStarList->front()->getCcdImage()) == &ccdImage);
 
   58     std::size_t nparModel = (_fittingModel) ? _photometryModel->getNpar(ccdImage) : 0;
 
   62     if (_fittingModel) _photometryModel->getMappingIndices(ccdImage, indices);
 
   64     Eigen::VectorXd H(nparTotal);  
 
   66     Eigen::Index kTriplets = tripletList.getNextFreeIndex();
 
   67     const MeasuredStarList &catalog = (measuredStarList) ? *measuredStarList : ccdImage.getCatalogForFit();
 
   69     for (
auto const &measuredStar : catalog) {
 
   70         if (!measuredStar->isValid()) 
continue;
 
   73         double residual = _photometryModel->computeResidual(ccdImage, *measuredStar);
 
   74         double inverseSigma = 1.0 / _photometryModel->transformError(ccdImage, *measuredStar);
 
   75         double W = 
std::pow(inverseSigma, 2);
 
   78             _photometryModel->computeParameterDerivatives(*measuredStar, ccdImage, H);
 
   80                 Eigen::Index l = indices[k];
 
   81                 tripletList.addTriplet(l, kTriplets, H[k] * inverseSigma);
 
   82                 grad[l] += H[k] * W * residual;
 
   86             Eigen::Index index = measuredStar->getFittedStar()->getIndexInMatrix();
 
   88             tripletList.addTriplet(index, kTriplets, -1.0 * inverseSigma);
 
   89             grad[index] += -1.0 * W * residual;
 
   94     tripletList.setNextFreeIndex(kTriplets);
 
   97 void PhotometryFit::leastSquareDerivativesReference(FittedStarList 
const &fittedStarList,
 
   98                                                     TripletList &tripletList, Eigen::VectorXd &grad)
 const {
 
  105     if (!_fittingFluxes) 
return;
 
  109     Eigen::Index kTriplets = tripletList.getNextFreeIndex();
 
  111     for (
auto const &fittedStar : fittedStarList) {
 
  112         auto refStar = fittedStar->getRefStar();
 
  113         if (refStar == 
nullptr) 
continue;  
 
  121         double inverseSigma = 1.0 / _photometryModel->getRefError(*refStar);
 
  123         double residual = _photometryModel->computeRefResidual(*fittedStar, *refStar);
 
  125         Eigen::Index index = fittedStar->getIndexInMatrix();
 
  127         tripletList.addTriplet(index, kTriplets, 1.0 * inverseSigma);
 
  128         grad(index) += 1.0 * 
std::pow(inverseSigma, 2) * residual;
 
  131     tripletList.setNextFreeIndex(kTriplets);
 
  134 void PhotometryFit::accumulateStatImageList(
CcdImageList const &ccdImageList, Chi2Accumulator &accum)
 const {
 
  139     for (
auto const &ccdImage : ccdImageList) {
 
  140         auto &catalog = ccdImage->getCatalogForFit();
 
  142         for (
auto const &measuredStar : catalog) {
 
  143             if (!measuredStar->isValid()) 
continue;
 
  144             double sigma = _photometryModel->transformError(*ccdImage, *measuredStar);
 
  145             double residual = _photometryModel->computeResidual(*ccdImage, *measuredStar);
 
  148             accum.addEntry(chi2Val, 1, measuredStar);
 
  153 void PhotometryFit::accumulateStatRefStars(Chi2Accumulator &accum)
 const {
 
  159     FittedStarList &fittedStarList = 
_associations->fittedStarList;
 
  160     for (
auto const &fittedStar : fittedStarList) {
 
  161         auto refStar = fittedStar->getRefStar();
 
  162         if (refStar == 
nullptr) 
continue;
 
  163         double sigma = _photometryModel->getRefError(*refStar);
 
  164         double residual = _photometryModel->computeRefResidual(*fittedStar, *refStar);
 
  166         accum.addEntry(chi2, 1, fittedStar);
 
  173 void PhotometryFit::getIndicesOfMeasuredStar(MeasuredStar 
const &measuredStar, 
IndexVector &indices)
 const {
 
  176         _photometryModel->getMappingIndices(measuredStar.getCcdImage(), indices);
 
  178     if (_fittingFluxes) {
 
  180         Eigen::Index fsIndex = fs->getIndexInMatrix();
 
  189     _fittingFluxes = (
_whatToFit.
find(
"Fluxes") != std::string::npos);
 
  192     _nModelParams = (_fittingModel) ? _photometryModel->assignIndices(whatToFit, 0) : 0;
 
  195     if (_fittingFluxes) {
 
  201             fittedStar->setIndexInMatrix(ipar);
 
  214                           "PhotometryFit::offsetParams : the provided vector length is not compatible with " 
  215                           "the current whatToFit setting");
 
  216     if (_fittingModel) _photometryModel->offsetParams(delta);
 
  218     if (_fittingFluxes) {
 
  223             Eigen::Index index = fittedStar->getIndexInMatrix();
 
  224             _photometryModel->offsetFittedStar(*fittedStar, delta(index));
 
  233     ofile << 
"id" << separator << 
"xccd" << separator << 
"yccd" << separator;
 
  234     ofile << 
"mag" << separator << 
"instMag" << separator << 
"instMagErr" << separator;
 
  235     ofile << 
"instFlux" << separator << 
"instFluxErr" << separator;
 
  236     ofile << 
"inputFlux" << separator << 
"inputFluxErr" << separator;
 
  237     ofile << 
"transformedFlux" << separator << 
"transformedFluxErr" << separator;
 
  238     ofile << 
"fittedFlux" << separator;
 
  239     ofile << 
"epoch" << separator;
 
  240     ofile << 
"fsindex" << separator;
 
  241     ofile << 
"ra" << separator << 
"dec" << separator;
 
  242     ofile << 
"chi2" << separator << 
"nm" << separator;
 
  243     ofile << 
"detector" << separator << 
"visit" << separator << 
std::endl;
 
  245     ofile << 
"id in source catalog" << separator << 
"coordinates in CCD" << separator << separator;
 
  246     ofile << 
"fitted magnitude" << separator << 
"measured magnitude" << separator
 
  247           << 
"measured magnitude error" << separator;
 
  248     ofile << 
"measured instrumental flux (ADU)" << separator << 
"measured instrument flux error" << separator;
 
  249     ofile << 
"measured flux (nJy)" << separator << 
"measured flux error" << separator;
 
  250     ofile << separator << separator;
 
  251     ofile << 
"fitted flux (nJy)" << separator;
 
  252     ofile << 
"Julian Epoch Year of the measurement" << separator;
 
  253     ofile << 
"unique index of the fittedStar" << separator;
 
  254     ofile << 
"on-sky position of fitted star" << separator << separator;
 
  255     ofile << 
"contribution to chi2 (1 dof)" << separator << 
"number of measurements of this FittedStar" 
  257     ofile << 
"detector id" << separator << 
"visit id" << 
std::endl;
 
  260     for (
auto const &ccdImage : ccdImageList) {
 
  262         for (
auto const &measuredStar : cat) {
 
  263             if (!measuredStar->isValid()) 
continue;
 
  265             double instFluxErr = _photometryModel->tweakFluxError(*measuredStar);
 
  266             double flux = _photometryModel->transform(*ccdImage, *measuredStar);
 
  267             double fluxErr = _photometryModel->transformError(*ccdImage, *measuredStar);
 
  269             double residual = _photometryModel->computeResidual(*ccdImage, *measuredStar);
 
  270             double chi2Val = 
std::pow(residual / fluxErr, 2);
 
  273             ofile << measuredStar->getId() << separator << measuredStar->x << separator << measuredStar->y
 
  275             ofile << fittedStar->getMag() << separator << measuredStar->getInstMag() << separator
 
  276                   << measuredStar->getInstMagErr() << separator;
 
  277             ofile << measuredStar->getInstFlux() << separator << instFluxErr << separator;
 
  278             ofile << measuredStar->getFlux() << separator << measuredStar->getFluxErr() << separator;
 
  279             ofile << flux << separator << fluxErr << separator << fittedStar->getFlux() << separator;
 
  280             ofile << ccdImage->getEpoch() << separator;
 
  281             ofile << fittedStar->getIndexInMatrix() << separator;
 
  282             ofile << fittedStar->x << separator << fittedStar->y << separator;
 
  283             ofile << chi2Val << separator << fittedStar->getMeasurementCount() << separator;
 
  284             ofile << ccdImage->getCcdId() << separator << ccdImage->getVisit() << 
std::endl;
 
  293     ofile << 
"ra" << separator << 
"dec " << separator;
 
  294     ofile << 
"mag" << separator;
 
  295     ofile << 
"refFlux" << separator << 
"refFluxErr" << separator;
 
  296     ofile << 
"fittedFlux" << separator << 
"fittedFluxErr" << separator;
 
  297     ofile << 
"fsindex" << separator << 
"chi2" << separator << 
"nm" << 
std::endl;
 
  299     ofile << 
"coordinates of fittedStar" << separator << separator;
 
  300     ofile << 
"magnitude" << separator;
 
  301     ofile << 
"refStar flux (nJy)" << separator << 
"refStar fluxErr" << separator;
 
  302     ofile << 
"fittedStar flux (nJy)" << separator << 
"fittedStar fluxErr" << separator;
 
  303     ofile << 
"unique index of the fittedStar" << separator << 
"refStar contribution to chi2 (1 dof)" 
  304           << separator << 
"number of measurements of this FittedStar" << 
std::endl;
 
  308     for (
auto const &fittedStar : fittedStarList) {
 
  309         const RefStar *refStar = fittedStar->getRefStar();
 
  310         if (refStar == 
nullptr) 
continue;
 
  315         ofile << fittedStar->x << separator << fittedStar->y << separator;
 
  316         ofile << fittedStar->getMag() << separator;
 
  318         ofile << fittedStar->getFlux() << separator << fittedStar->getFluxErr() << separator;
 
  319         ofile << fittedStar->getIndexInMatrix() << separator << chi2 << separator
 
  320               << fittedStar->getMeasurementCount() << 
std::endl;
 
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
 
afw::table::Key< double > sigma
 
LSST DM logging module built on log4cxx.
 
#define LOGLS_INFO(logger, message)
Log a info-level message using an iostream-based interface.
 
#define LOGLS_DEBUG(logger, message)
Log a debug-level message using an iostream-based interface.
 
double getFluxErr() const
 
A list of FittedStar s. Such a list is typically constructed by Associations.
 
Eigen::Index _nModelParams
 
std::shared_ptr< Associations > _associations
 
Eigen::Index _nStarParams
 
A list of MeasuredStar. They are usually filled in Associations::createCcdImage.
 
void assignIndices(std::string const &whatToFit) override
Set parameters to fit and assign indices in the big matrix.
 
void saveChi2RefContributions(std::string const &filename) const override
Save a CSV file containing residuals of reference terms.
 
void saveChi2MeasContributions(std::string const &filename) const override
Save a CSV file containing residuals of measurement terms.
 
void offsetParams(Eigen::VectorXd const &delta) override
Offset the parameters by the requested quantities.
 
Objects used as position/flux anchors (e.g.
 
Reports invalid arguments.
 
std::list< std::shared_ptr< CcdImage > > CcdImageList
 
A base class for image defects.
 
T setprecision(T... args)