9#ifndef __KLIPreduction_hpp__
10#define __KLIPreduction_hpp__
17#include "../ipc/ompLoopWatcher.hpp"
18#include "../math/eigenLapack.hpp"
19#include "../math/geo.hpp"
20#include "../sigproc/gramSchmidt.hpp"
21using namespace mx::sigproc;
49std::string excludeMethodStr(
int method );
51int excludeMethodFmStr(
const std::string &method );
65std::string includeMethodStr(
int method );
67int includeMethodFmStr(
const std::string &method );
82template <
typename _realT,
class _derotFunctObj,
typename _evCalcT =
double>
87 typedef Eigen::Array<realT, Eigen::Dynamic, Eigen::Dynamic> eigenImageT;
89 typedef _evCalcT evCalcT;
100 const std::string &prefix,
101 const std::string &ext
119 const std::string &prefix,
120 const std::string &ext,
121 const std::string &RDIdir,
122 const std::string &RDIprefix,
124 const std::string &RDIext =
""
139 const std::string &RDIfileListFile
186 int m_maxNmodes{ 0 };
246 bool m_rightReason =
false;
248 realT m_rightReasonRadius = 2.5;
259 std::vector<realT> &sds
263 std::vector<_realT> m_minr;
264 std::vector<_realT> m_maxr;
265 std::vector<_realT> m_minq;
266 std::vector<_realT> m_maxq;
276 int regions(
const std::vector<realT> & minr,
277 const std::vector<realT> & maxr,
278 const std::vector<realT> & minq,
279 const std::vector<realT> & maxq
295 std::vector<realT> vminr( 1, minr );
296 std::vector<realT> vmaxr( 1, maxr );
297 std::vector<realT> vminq( 1, minq );
298 std::vector<realT> vmaxq( 1, maxq );
300 return regions( vminr, vmaxr, vminq, vmaxq );
306 std::vector<size_t> &idx,
312 int processPSFSub(
const std::string &dir,
const std::string &prefix,
const std::string &ext );
314 double t_worker_begin{ 0 };
315 double t_worker_end{ 0 };
317 double t_eigenv{ 0 };
323 printf(
"KLIP reduction times: \n" );
324 printf(
" Total time: %f sec\n", this->t_end - this->t_begin );
325 printf(
" Loading: %f sec\n", this->t_load_end - this->t_load_begin );
326 printf(
" Fake Injection: %f sec\n", this->t_fake_end - this->t_fake_begin );
327 printf(
" Coadding: %f sec\n", this->t_coadd_end - this->t_coadd_begin );
328 printf(
" Preprocessing: %f sec\n", this->t_preproc_end - this->t_preproc_begin );
329 printf(
" Az USM: %f sec\n", this->t_azusm_end - this->t_azusm_begin );
330 printf(
" Gauss USM: %f sec\n", this->t_gaussusm_end - this->t_gaussusm_begin );
331 printf(
" KLIP algorithm: %f elapsed real sec\n", this->t_worker_end - this->t_worker_begin );
332 double klip_cpu = this->t_eigenv + this->t_klim + this->t_psf;
333 printf(
" EigenDecomposition %f cpu sec (%f%%)\n", this->t_eigenv, this->t_eigenv / klip_cpu * 100 );
334 printf(
" KL image calc %f cpu sec (%f%%)\n", this->t_klim, this->t_klim / klip_cpu * 100 );
335 printf(
" PSF calc/sub %f cpu sec (%f%%)\n", this->t_psf, this->t_psf / klip_cpu * 100 );
336 printf(
" Derotation: %f sec\n", this->t_derotate_end - this->t_derotate_begin );
337 printf(
" Combination: %f sec\n", this->t_combo_end - this->t_combo_begin );
341template <
typename _realT,
class _derotFunctObj,
typename _evCalcT>
346template <
typename _realT,
class _derotFunctObj,
typename _evCalcT>
348 const std::string &prefix,
349 const std::string &ext )
354template <
typename _realT,
class _derotFunctObj,
typename _evCalcT>
360template <
typename _realT,
class _derotFunctObj,
typename _evCalcT>
362 const std::string &prefix,
363 const std::string &ext,
364 const std::string &RDIdir,
365 const std::string &RDIprefix,
366 const std::string &RDIext )
367 :
ADIobservation<_realT, _derotFunctObj>( dir, prefix, ext, RDIdir, RDIprefix, RDIext )
371template <
typename _realT,
class _derotFunctObj,
typename _evCalcT>
373 const std::string &RDIfileListFile )
374 :
ADIobservation<_realT, _derotFunctObj>( fileListFile, RDIfileListFile )
378template <
typename _realT,
class _derotFunctObj,
typename _evCalcT>
383template <
typename realT,
class derotFunctObj,
typename evCalcT>
387 std::vector<realT> &norms )
390 norms.resize( rims.planes() );
392 bool haveMask =
false;
393 if( cmask.rows() > 0 && cmask.cols() > 0 )
411 for(
int n = 0; n < rims.planes(); ++n )
413 rims.
image( n ) -= mean;
417 rims.
image( n ) *= cmask;
420 realT immean = rims.
image( n ).mean();
421 norms[n] = ( rims.
image( n ) - immean ).matrix().norm();
435 for(
int n = 0; n < tims.planes(); ++n )
437 tims.
image( n ) -= mean;
441 tims.
image( n ) *= cmask;
449 std::vector<realT> work;
451 for(
int n = 0; n < rims.planes(); ++n )
455 mean = rims.
image( n ).mean();
462 rims.
image( n ) -= mean;
466 rims.
image( n ) *= cmask;
470 realT immean = rims.
image( n ).mean();
471 norms[n] = ( rims.
image( n ) - immean ).matrix().norm();
476 for(
int n = 0; n < tims.planes(); ++n )
480 mean = tims.
image( n ).mean();
487 tims.
image( n ) -= mean;
490 tims.
image( n ) *= cmask;
498 mxThrowException(
err::invalidconfig,
"KlipReduction::meanSubtract",
"pixelTSNormMethod is unknown");
502 mxThrowException(
err::invalidconfig,
"KlipReduction::meanSubtract",
"pixelTSNormMethod is rmsSigmaClipped, which is not implemented");
507 std::cerr <<
"normalizing pixels\n";
508 std::vector<realT> pixs(rims.planes());
510 for(
int cc = 0; cc < rims.cols(); ++cc)
512 for(
int rr = 0; rr < rims.rows(); ++rr)
516 if(cmask(rr,cc) == 0)
523 for(
int pp = 0; pp < rims.planes(); ++pp)
525 pixs[pp] = rims.
image(pp)(rr,cc);
530 for(
int pp = 0; pp < rims.planes(); ++pp)
532 rims.
image(pp)(rr,cc) /= sd;
539template <
typename _realT,
class _derotFunctObj,
typename _evCalcT>
541 const std::vector<_realT> & maxr,
542 const std::vector<_realT> & minq,
543 const std::vector<_realT> & maxq )
552 m_maxNmodes = m_Nmodes[0];
553 for(
size_t i = 1; i < m_Nmodes.size(); ++i )
555 if( m_Nmodes[i] > m_maxNmodes )
556 m_maxNmodes = m_Nmodes[i];
559 std::cerr <<
"Beginning\n";
561 if( this->m_imSize == 0 )
563 this->m_imSize = 2 * ( *std::max_element( m_maxr.begin(), m_maxr.end() ) + m_padSize );
565 std::cerr <<
"set image size based on regions to " << this->m_imSize <<
"\n";
568 if( !this->m_filesRead )
570 if( this->readFiles() < 0 )
575 if( !this->m_RDIfilesRead && this->m_RDIfileList.size() != 0 )
577 if( this->readRDIFiles() < 0 )
583 if( this->m_preProcess_only && !this->m_skipPreProcess )
585 std::cerr <<
"Pre-processing complete, stopping.\n";
589 std::cerr <<
"allocating psf subtracted cubes\n";
590 this->m_psfsub.resize( m_Nmodes.size() );
591 for(
size_t n = 0; n < m_Nmodes.size(); ++n )
593 this->m_psfsub[n].resize( this->m_Nrows, this->m_Ncols, this->m_Nims );
594 this->m_psfsub[n].cube().setZero();
598 eigenImageT rIm( this->m_Nrows, this->m_Ncols );
599 eigenImageT qIm( this->m_Nrows, this->m_Ncols );
601 radAngImage<math::degreesT<realT>>( rIm, qIm, .5 * ( this->m_Nrows - 1 ), .5 * ( this->m_Ncols - 1 ) );
603 m_imsIncluded.resize( this->m_Nims, this->m_Nims );
604 m_imsIncluded.setConstant( 1 );
606 if( this->m_refIms.planes() > 0 )
608 std::cerr <<
"******* RDI MODE **********\n";
612 std::cerr <<
"******* ADI MODE **********\n";
615 std::cerr <<
"processing " << minr.size() <<
" regions\n";
618 for(
size_t regno = 0; regno < minr.size(); ++regno )
620 std::cerr <<
" region " << regno+1 <<
": " << m_minr[regno] <<
"-" << m_maxr[regno] <<
" pixels, ";
621 std::cerr << m_minq[regno] <<
"-" << m_maxq[regno] <<
" degrees. \n";
623 eigenImageT *maskPtr =
nullptr;
625 if( this->m_mask.rows() == this->m_Nrows && this->m_mask.cols() == this->m_Ncols )
627 maskPtr = &this->m_mask;
630 std::vector<size_t> idx = annulusIndices<math::degreesT<realT>>( rIm,
632 .5 * ( this->m_Nrows - 1 ),
633 .5 * ( this->m_Ncols - 1 ),
649 if( this->m_refIms.planes() > 0 )
651 rims.resize( idx.size(), 1, this->m_Nims );
654 for(
int i = 0; i < this->m_Nims; ++i )
656 auto tim = tims.image( i );
660 for(
int p = 0; p < this->m_refIms.planes(); ++p )
662 auto rim = rims.
image( p );
666 if( this->m_maskFile !=
"" )
673 std::vector<realT> sds;
676 meanSubtract( tims, tims, cmask, sds );
680 ffff.
write(
"tims.fits", tims);
682 std::ofstream fout(
"idx.dat");
683 std::ofstream aout(
"derot.dat");
690 for(
int pp =0; pp < tims.planes(); ++pp)
692 aout << this->m_derotF.derotAngle( pp ) <<
"\n";
696 fout.open(
"dims.dat");
697 fout << this->m_Nrows <<
"\n";
698 fout << this->m_Ncols <<
"\n";
718 if( this->m_refIms.planes() > 0 )
726 dang = fabs( atan( m_minDPx / minr[regno] ) );
739 dangMax = fabs( atan( m_maxDPx / minr[regno] ) );
752 if( this->m_refIms.planes() > 0 )
754 worker( rims, tims, cmask, idx, dang, dangMax );
758 worker( tims, tims, cmask, idx, dang, dangMax );
763 ffii.
write(
"imsIncluded.fits", m_imsIncluded );
765 if( finalProcess() < 0 )
767 std::cerr <<
"Error in final processing\n";
782 bool included{
true };
785template <
typename eigenT,
typename eigenTin>
786void extractRowsAndCols( eigenT &out,
const eigenTin &in,
const std::vector<size_t> &idx )
788 out.resize( idx.size(), idx.size() );
790 for(
size_t i = 0; i < idx.size(); ++i )
792 for(
size_t j = 0; j < idx.size(); ++j )
794 out( i, j ) = in( idx[i], idx[j] );
799template <
typename eigenT,
typename eigenTin>
800void extractCols( eigenT &out,
const eigenTin &in,
const std::vector<size_t> &idx )
802 out.resize( in.rows(), idx.size() );
804 for(
size_t i = 0; i < idx.size(); ++i )
806 out.col( i ) = in.col( idx[i] );
810template <
typename realT,
typename eigenT,
typename eigenTv,
class derotFunctObj>
811void collapseCovar( eigenT &cutCV,
813 const std::vector<realT> &sds,
821 int excludeMethodMax,
823 const derotFunctObj &derotF,
824 eigenImage<int> &imsIncluded )
826 std::vector<cvEntry> allidx( Nims );
829 for(
int i = 0; i < Nims; ++i )
832 allidx[i].angle = derotF.derotAngle( i );
837 allidx[i].cvVal = CV( imno, i ) / ( sds[i] * sds[imno] );
841 allidx[i].cvVal = CV( i, imno ) / ( sds[i] * sds[imno] );
847 for(
size_t j = 0; j < Nims; ++j )
852 allidx[j].included =
false;
858 for(
size_t j = 0; j < Nims; ++j )
860 if( fabs( (
long)j - imno ) <= dang )
862 allidx[j].included =
false;
869 for(
size_t j = 0; j < Nims; ++j )
873 allidx[j].included =
false;
878 for(
size_t j = 0; j < Nims; ++j )
880 if( fabs( (
long)j - imno ) > dangMax )
881 allidx[j].included =
false;
885 if( includeRefNum > 0 && (
size_t)includeRefNum < allidx.size() )
888 for(
size_t j = 0; j < Nims; ++j )
890 if( allidx[j].included ==
true )
897 std::vector<realT> cvVal;
898 cvVal.resize( kept );
900 for(
size_t j = 0; j < Nims; ++j )
902 if( allidx[j].included ==
true )
904 cvVal[
k] = allidx[j].cvVal;
910 std::nth_element( cvVal.begin(), cvVal.begin() + ( kept - includeRefNum ), cvVal.end() );
912 realT mincorr = cvVal[kept - includeRefNum];
915 for(
size_t j = 0; j < Nims; ++j )
917 if( allidx[j].cvVal < mincorr )
919 allidx[j].included =
false;
924 std::vector<size_t> keepidx;
925 for(
size_t j = 0; j < Nims; ++j )
927 imsIncluded( imno, j ) = allidx[j].included;
929 if( allidx[j].included )
930 keepidx.push_back( j );
937 if( keepidx.size() == 0 )
939 std::cerr <<
"\n\n" << imno <<
"\n\n";
942 extractRowsAndCols( cutCV, CV, keepidx );
943 extractCols( rimsCut, rims, keepidx );
946template <
typename _realT,
class _derotFunctObj,
typename _evCalcT>
950 std::vector<size_t> &idx,
957 std::vector<realT> sds;
962 meanSubtract( rims, tims, cmask, sds );
970 ff.
write(
"cv.fits", cv );
974 eigenImageT master_klims;
975 eigenImageT master_projMat;
980 rrMask.resize( tims.rows(), tims.rows() );
981 rrMask.setConstant( 1 );
984 rrmim.resize( this->m_Nrows, this->m_Ncols );
985 for(
int rr = 0; rr < rrMask.rows(); ++rr )
987 rrmim.setConstant( 1 );
990 int jj = idx[rr] / this->m_Ncols;
991 int ii = idx[rr] - this->m_Ncols * jj;
993 maskCircle( rrmim, ii, jj, m_rightReasonRadius, 0, 0 );
996 for(
int cc = 0; cc < rrMask.cols(); ++cc )
998 rrMask( rr, cc ) = rrmim.data()[idx[cc]];
1002 ff.
write(
"rrMask.fits", rrMask );
1014 master_projMat = ( master_klims.matrix().transpose() * master_klims.matrix() ).array();
1015 ff.
write(
"projMat.fits", master_projMat );
1016 master_projMat *= rrMask;
1017 ff.
write(
"projMatrr.fits", master_projMat );
1020 t_eigenv += teigenv;
1025 #pragma omp parallel
1032 eigenImageT rims_cut;
1035 eigenImageT projMat;
1040 m_includeRefNum == 0 )
1042 klims = master_klims;
1045 projMat = master_projMat;
1052 for(
int imno = 0; imno < this->m_Nims; ++imno )
1056 collapseCovar<realT>( cv_cut,
1065 this->m_excludeMethod,
1066 this->m_excludeMethodMax,
1067 this->m_includeRefNum,
1072 double teigenv, tklim;
1073 math::calcKLModes( klims, cv_cut, rims_cut, m_maxNmodes, &mem, &teigenv, &tklim );
1077 projMat = ( klims.matrix().transpose() * klims.matrix() ).array();
1081 t_eigenv += teigenv;
1084 cfs.resize( 1, klims.rows() );
1088 if( !m_rightReason )
1090 for(
int j = 0; j < cfs.size(); ++j )
1092 cfs( j ) = klims.row( j ).matrix().dot( tims.
cube().col( imno ).matrix() );
1095 for(
size_t mode_i = 0; mode_i < m_Nmodes.size(); ++mode_i )
1097 psf = cfs( cfs.size() - 1 ) * klims.row( cfs.size() - 1 );
1103 for(
int j = cfs.size() - 2; j >= cfs.size() - m_Nmodes[mode_i] && j >= 0; --j )
1105 psf += cfs( j ) * klims.row( j );
1110 this->m_psfsub[mode_i].cube().col( imno ), tims.
cube().col( imno ) - psf.transpose(), idx );
1115 psf = projMat.matrix() * tims.
cube().col( imno ).matrix();
1129template <
typename _realT,
class _derotFunctObj,
typename _evCalcT>
1132 if( this->m_postMedSub )
1134 std::cerr <<
"Subtracting medians in post\n";
1136 for(
size_t n = 0; n < this->m_psfsub.size(); ++n )
1139 #pragma omp parallel
1144 this->m_psfsub[n].median( medim );
1149 for(
int i = 0; i < this->m_psfsub[n].planes(); ++i )
1151 this->m_psfsub[n].image( i ) -= medim;
1157 if( this->m_doDerotate )
1159 std::cerr <<
"derotating\n";
1163 if( this->m_combineMethod > 0 )
1165 std::cerr <<
"combining\n";
1166 this->combineFinim();
1169 if( this->m_doWriteFinim ==
true || this->m_doOutputPSFSub ==
true )
1171 std::cerr <<
"writing\n";
1175 this->ADIobservation<_realT, _derotFunctObj>::stdFitsHeader( &head );
1177 head.
append(
"", fits::fitsCommentType(),
"----------------------------------------" );
1178 head.
append(
"", fits::fitsCommentType(),
"mx::KLIPreduction parameters:" );
1179 head.
append(
"", fits::fitsCommentType(),
"----------------------------------------" );
1181 head.
append(
"MEAN SUB METHOD", HCI::meanSubMethodStr( m_meanSubMethod ),
"PCA mean subtraction method" );
1182 head.
append(
"PIXTS NORM METHOD", HCI::pixelTSNormMethodStr( m_pixelTSNormMethod ),
"Pixel TS norm method" );
1184 std::stringstream str;
1186 if( m_Nmodes.size() > 0 )
1188 for(
size_t nm = 0; nm < m_Nmodes.size() - 1; ++nm )
1189 str << m_Nmodes[nm] <<
",";
1190 str << m_Nmodes[m_Nmodes.size() - 1];
1191 head.
append<
char *>(
"NMODES", (
char *)str.str().c_str(),
"number of modes" );
1194 head.
append<
bool>(
"RIGHT REASON", m_rightReason,
"whether or not the right reason mask is used");
1197 head.
append<realT>(
"RIGHT REASON RADIUS", m_rightReasonRadius,
"radius of the right reason mask");
1200 if( m_minr.size() > 0 )
1203 for(
size_t nm = 0; nm < m_minr.size() - 1; ++nm )
1204 str << m_minr[nm] <<
",";
1205 str << m_minr[m_minr.size() - 1];
1206 head.
append<
char *>(
"REGMINR", (
char *)str.str().c_str(),
"region inner edge(s)" );
1209 if( m_maxr.size() > 0 )
1212 for(
size_t nm = 0; nm < m_maxr.size() - 1; ++nm )
1213 str << m_maxr[nm] <<
",";
1214 str << m_maxr[m_maxr.size() - 1];
1215 head.
append<
char *>(
"REGMAXR", (
char *)str.str().c_str(),
"region outer edge(s)" );
1218 if( m_minq.size() > 0 )
1221 for(
size_t nm = 0; nm < m_minq.size() - 1; ++nm )
1222 str << m_minq[nm] <<
",";
1223 str << m_minq[m_minq.size() - 1];
1224 head.
append<
char *>(
"REGMINQ", (
char *)str.str().c_str(),
"region minimum angle(s)" );
1227 if( m_maxq.size() > 0 )
1230 for(
size_t nm = 0; nm < m_maxq.size() - 1; ++nm )
1231 str << m_maxq[nm] <<
",";
1232 str << m_maxq[m_maxq.size() - 1];
1233 head.
append<
char *>(
"REGMAXQ", (
char *)str.str().c_str(),
"region maximum angle(s)" );
1236 head.
append<std::string>(
"EXMTHDMN", HCI::excludeMethodStr( m_excludeMethod ),
"exclusion method (min)" );
1237 head.
append<realT>(
"MINDPX", m_minDPx,
"minimum delta (units based on EXMTHDMN)" );
1239 head.
append<std::string>(
"EXMTHDMX", HCI::excludeMethodStr( m_excludeMethodMax ),
"exclusion method (max)" );
1240 head.
append<realT>(
"MAXDPX", m_maxDPx,
"maximum delta (units based on EXMTHDMX)" );
1242 head.
append<std::string>(
"INMTHDMX", HCI::includeMethodStr( m_includeMethod ),
"inclusion method" );
1243 head.
append<
int>(
"INCLREFN", m_includeRefNum,
"number of images included by INMTHDMX" );
1245 if( this->m_doWriteFinim ==
true && this->m_combineMethod > 0 )
1247 this->writeFinim( &head );
1250 if( this->m_doOutputPSFSub )
1252 this->outputPSFSub( &head );
1259template <
typename _realT,
class _derotFunctObj,
typename _evCalcT>
1260int KLIPreduction<_realT, _derotFunctObj, _evCalcT>::processPSFSub(
const std::string &dir,
1261 const std::string &prefix,
1262 const std::string &ext )
1265 std::cerr <<
"Beginning PSF Subtracted Image Processing\n";
1271 eigenImage<realT> im;
1274 ff.
read( im, fh, flist[0] );
1276 if( fh.
count(
"MEANSUBM" ) == 0 )
1278 mxError(
"KLIPReduction", MXE_PARAMNOTSET,
"MEANSUBM not found in FITS header." );
1281 m_meanSubMethod = HCI::meanSubMethodFmStr( fh[
"MEANSUBM"].String() );
1282 std::cerr <<
"meanSubMethod: " << HCI::meanSubMethodStr( m_meanSubMethod ) <<
"\n";
1284 if( fh.
count(
"NMODES" ) == 0 )
1286 mxError(
"KLIPReduction", MXE_PARAMNOTSET,
"NMODES not found in FITS header." );
1290 if( m_Nmodes.size() == 0 )
1292 mxError(
"KLIPReduction", MXE_PARSEERR,
"NMODES vector did not parse correctly." );
1295 std::cerr <<
"nModes: " << fh[
"NMODES"].String() <<
"\n";
1298 if( fh.
count(
"REGMINR" ) == 0 )
1300 mxError(
"KLIPReduction", MXE_PARAMNOTSET,
"REGMINR not found in FITS header." );
1304 if( m_minr.size() == 0 )
1306 mxError(
"KLIPReduction", MXE_PARSEERR,
"REGMINR vector did not parse correctly." );
1309 std::cerr <<
"minr: " << fh[
"REGMINR"].String() <<
"\n";
1312 if( fh.
count(
"REGMAXR" ) == 0 )
1314 mxError(
"KLIPReduction", MXE_PARAMNOTSET,
"REGMAXR not found in FITS header." );
1318 if( m_maxr.size() == 0 )
1320 mxError(
"KLIPReduction", MXE_PARSEERR,
"REGMAXR vector did not parse correctly." );
1323 std::cerr <<
"minr: " << fh[
"REGMAXR"].String() <<
"\n";
1326 if( fh.
count(
"REGMINQ" ) == 0 )
1328 mxError(
"KLIPReduction", MXE_PARAMNOTSET,
"REGMINQ not found in FITS header." );
1332 if( m_minq.size() == 0 )
1334 mxError(
"KLIPReduction", MXE_PARSEERR,
"REGMINQ vector did not parse correctly." );
1337 std::cerr <<
"minq: " << fh[
"REGMINQ"].String() <<
"\n";
1340 if( fh.
count(
"REGMAXR" ) == 0 )
1342 mxError(
"KLIPReduction", MXE_PARAMNOTSET,
"REGMAXR not found in FITS header." );
1346 if( m_maxq.size() == 0 )
1348 mxError(
"KLIPReduction", MXE_PARSEERR,
"REGMAXR vector did not parse correctly." );
1351 std::cerr <<
"minr: " << fh[
"REGMAXR"].String() <<
"\n";
1353 if( fh.
count(
"EXMTHDMN" ) == 0 )
1355 mxError(
"KLIPReduction", MXE_PARAMNOTSET,
"EXMTHDMN not found in FITS header." );
1358 m_excludeMethod = HCI::excludeMethodFmStr( fh[
"EXMTHDMN"].String() );
1359 std::cerr <<
"excludeMethod: " << HCI::excludeMethodStr( m_excludeMethod ) <<
"\n";
1361 if( fh.
count(
"MINDPX" ) == 0 )
1363 mxError(
"KLIPReduction", MXE_PARAMNOTSET,
"MINDPX not found in FITS header." );
1366 m_minDPx = fh[
"MINDPX"].value<realT>();
1367 std::cerr <<
"minDPx: " << m_minDPx <<
"\n";
1369 if( fh.
count(
"EXMTHDMX" ) == 0 )
1371 mxError(
"KLIPReduction", MXE_PARAMNOTSET,
"EXMTHDMX not found in FITS header." );
1374 m_excludeMethodMax = HCI::excludeMethodFmStr( fh[
"EXMTHDMX"].String() );
1375 std::cerr <<
"excludeMethodMax: " << HCI::excludeMethodStr( m_excludeMethodMax ) <<
"\n";
1377 if( fh.
count(
"MAXDPX" ) == 0 )
1379 mxError(
"KLIPReduction", MXE_PARAMNOTSET,
"MAXDPX not found in FITS header." );
1382 m_maxDPx = fh[
"MAXDPX"].value<realT>();
1383 std::cerr <<
"maxDPx: " << m_maxDPx <<
"\n";
1385 if( fh.
count(
"INMTHDMX" ) == 0 )
1387 mxError(
"KLIPReduction", MXE_PARAMNOTSET,
"INMTHDMX not found in FITS header." );
1390 m_includeMethod = HCI::includeMethodFmStr( fh[
"INMTHDMX"].String() );
1391 std::cerr <<
"includeMethod: " << HCI::includeMethodStr( m_includeMethod ) <<
"\n";
1393 if( fh.
count(
"INCLREFN" ) == 0 )
1395 mxError(
"KLIPReduction", MXE_PARAMNOTSET,
"INCLREFN not found in FITS header." );
1398 m_includeRefNum = fh[
"INCLREFN"].value<
int>();
1399 std::cerr <<
"includedRefNum: " << m_includeRefNum <<
"\n";
1401 this->m_skipPreProcess =
true;
1403 this->m_keywords.
clear();
1405 this->readPSFSub( dir, prefix, ext, m_Nmodes.size() );
1414template <
typename realT>
1417extern template struct KLIPreduction<float, ADIDerotator<float>,
float>;
1418extern template struct KLIPreduction<float, ADIDerotator<float>,
double>;
1419extern template struct KLIPreduction<double, ADIDerotator<double>,
double>;
Defines the ADI high contrast imaging data type.
mxException for invalid config settings
Class to manage interactions with a FITS file.
int write(const dataT *im, int d1, int d2, int d3, fitsHeader *head)
Write the contents of a raw array to the FITS file.
int read(dataT *data)
Read the contents of the FITS file into an array.
An image cube with an Eigen-like API.
void median(eigenT &mim)
Calculate the median image of the cube.
Eigen::Map< Eigen::Array< dataT, Eigen::Dynamic, Eigen::Dynamic > > cube()
Returns a 2D Eigen::Eigen::Map pointed at the entire cube.
void mean(eigenT &mim)
Calculate the mean image of the cube.
Eigen::Map< Eigen::Array< dataT, Eigen::Dynamic, Eigen::Dynamic > > asVectors()
Return an Eigen::Eigen::Map of the cube where each image is a vector.
Eigen::Map< Eigen::Array< dataT, Eigen::Dynamic, Eigen::Dynamic > > image(Index n)
Returns a 2D Eigen::Eigen::Map pointed at the specified image.
A class to track the number of iterations in an OMP parallelized loop.
void incrementAndOutputStatus()
Increment and output status.
constexpr units::realT k()
Boltzmann Constant.
imageT::Scalar imageMedian(const imageT &mat, const maskT *mask, std::vector< typename imageT::Scalar > *work=0)
Calculate the median of an Eigen-like array.
Eigen::Array< scalarT, -1, -1 > eigenImage
Definition of the eigenImage type, which is an alias for Eigen::Array.
MXLAPACK_INT calcKLModes(eigenT &klModes, eigenT &cv, const eigenT1 &Rims, int n_modes=0, syevrMem< _evCalcT > *mem=0, double *t_eigenv=nullptr, double *t_klim=nullptr)
Calculate the K-L modes, or principle components, given a covariance matrix.
void eigenSYRK(eigenT1 &cv, const eigenT2 &ims)
Calculates the lower triangular part of the covariance matrix of ims.
std::vector< std::string > getFileNames(const std::string &directory, const std::string &prefix, const std::string &substr, const std::string &extension)
Get a list of file names from the specified directory, specifying a prefix, a substring to match,...
constexpr floatT six_fifths()
Return 6/5 in the specified precision.
angleT::realT angleDiff(typename angleT::realT q1, typename angleT::realT q2)
Calculate the difference between two angles, correctly across 0/360.
realT dtor(realT q)
Convert from degrees to radians.
includeMethods
Image inclusion methods.
excludeMethods
Image exclusion methods.
meanSubMethod
Mean subtraction methods.
@ includeImno
include images which are closest in imno to the target
@ includeTime
include images which are closest in time to the target
@ includeAngle
include images which are closest in angle to the target
@ includeAll
include all images
@ includeCorr
include images which are most correlated with the target
@ excludePixel
Exclude by pixels of rotation.
@ excludeAngle
Exclude by angle of roration.
@ excludeImno
Exclude by number of images.
@ excludeNone
Exclude no images.
@ medianImage
The median image of the data is subtracted from each image.
@ meanImage
The mean image of the data is subtracted from each image.
void maskCircle(arrayT &m, typename arrayT::Scalar xcen, typename arrayT::Scalar ycen, typename arrayT::Scalar rad, typename arrayT::Scalar val, typename arrayT::Scalar pixbuf=0.5)
Mask a circle in an image.
void cutImageRegion(imageTout &imout, const imageTin &imin, const std::vector< size_t > &idx, bool resize=true)
Cut out a region of an image specified by an index-mask.
void insertImageRegion(imageTout imout, const imageTin &imin, const std::vector< size_t > &idx)
Insert a region of an image specified by an index-mask.
void parseStringVector(std::vector< typeT > &v, const std::string &s, char delim=',')
Parses a string into a vector of tokens delimited by a character.
typeT get_curr_time()
Get the current system time in seconds.
valueT vectorVariance(const valueT *vec, size_t sz, valueT mean)
Calculate the variance of a vector relative to a supplied mean value.
@ rmsSigmaClipped
the sigma clipped rms of the pixel time series
@ none
no pixel time series norm
@ unknown
unknown value, an error
Process an angular differential imaging (ADI) observation.
An implementation of the Karhunen-Loeve Image Processing (KLIP) algorithm.
HCI::meanSubMethod m_meanSubMethod
Specify how the data are centered for PCA within each search region.
HCI::pixelTSNormMethod m_pixelTSNormMethod
Specify if each pixel time-series is normalized.
void worker(eigenCube< realT > &rims, eigenCube< realT > &tims, eigenImageT &cmask, std::vector< size_t > &idx, realT dang, realT dangMax)
void meanSubtract(eigenCube< realT > &rims, eigenCube< realT > &tims, eigenImageT &cmask, std::vector< realT > &sds)
Subtract the basis mean from each of the images.
int m_includeRefNum
Number of reference images to include in the covariance matrix.
int regions(const std::vector< realT > &minr, const std::vector< realT > &maxr, const std::vector< realT > &minq, const std::vector< realT > &maxq)
Run KLIP in a set of geometric search regions.
std::vector< int > m_Nmodes
Specifies the number of modes to include in the PSF.
realT m_pixelTSSigma
Sigma-clipping parameter for pixel time-series normalization.
int regions(realT minr, realT maxr, realT minq, realT maxq)
Run KLIP in a geometric search region.
KLIPreduction()
Default c'tor.
int m_includeMethod
Controls how number of included images is calculated.
Type holding constants related to angle calculations in degrees.
A struct to hold the working memory for eigenSYEVR and maintain it between calls if desired.