9#ifndef __HCIobservation_hpp__
10#define __HCIobservation_hpp__
19#include "../mxlib.hpp"
21#include "../mxException.hpp"
23#include "../math/templateBLAS.hpp"
24#include "../sys/timeUtils.hpp"
25#include "../ioutils/fileUtils.hpp"
26#include "../ioutils/readColumns.hpp"
27#include "../ioutils/fits/fitsFile.hpp"
28#include "../ipc/ompLoopWatcher.hpp"
117template <
typename _realT>
125 typedef Eigen::Array<realT, Eigen::Dynamic, Eigen::Dynamic>
eigenImageT;
139 const std::string &prefix,
140 const std::string &ext
157 const std::string &prefix,
158 const std::string &ext,
160 const std::string &RDIdir,
161 const std::string &RDIprefix,
163 const std::string &RDIext =
""
176 const std::string &RDIfileListFile
185 const std::string &prefix,
186 const std::string &ext
203 const std::string &prefix,
204 const std::string &ext
429 int threshold( std::vector<std::string> &fileList,
430 const std::string &qualityFile,
431 realT qualityThreshold
467 std::vector<std::string> &coaddKeywords,
468 std::vector<double> &imageMJD,
469 std::vector<fits::fitsHeader> &heads,
699 int readPSFSub(
const std::string &dir,
const std::string &prefix,
const std::string &ext,
size_t nReductions );
704 double t_load_begin{ 0 };
705 double t_load_end{ 0 };
707 double t_coadd_begin{ 0 };
708 double t_coadd_end{ 0 };
710 double t_preproc_begin{ 0 };
711 double t_preproc_end{ 0 };
713 double t_azusm_begin{ 0 };
714 double t_azusm_end{ 0 };
716 double t_gaussusm_begin{ 0 };
717 double t_gaussusm_end{ 0 };
719 double t_combo_begin{ 0 };
720 double t_combo_end{ 0 };
725template <
typename _realT>
730template <
typename _realT>
733 load_fileList( dir, prefix, ext );
736template <
typename _realT>
739 load_fileList( fileListFile );
742template <
typename _realT>
744 const std::string &prefix,
745 const std::string &ext,
746 const std::string &RDIdir,
747 const std::string &RDIprefix,
748 const std::string &RDIext )
750 std::cerr <<
"HCI 6\n";
752 load_fileList( dir, prefix, ext );
760 load_RDIfileList( RDIdir, RDIprefix, re );
763template <
typename _realT>
766 load_fileList( fileListFile );
767 load_RDIfileList( RDIfileListFile );
770template <
typename _realT>
775 m_filesDeleted =
false;
778template <
typename _realT>
782 m_filesDeleted =
false;
785template <
typename _realT>
787 const std::string &prefix,
788 const std::string &ext )
791 m_RDIfilesDeleted =
false;
794template <
typename _realT>
798 m_RDIfilesDeleted =
false;
803template <
typename realT>
806 if( m_fileList.size() == 0 )
809 "HCIobservation<realT>::readFiles",
810 "The target fileList has 0 length, there are no files to be read." );
814 if( !m_filesDeleted )
816 if( m_deleteFront > 0 )
818 m_fileList.erase( m_fileList.begin(), m_fileList.begin() + m_deleteFront );
821 if( m_deleteBack > 0 )
823 m_fileList.erase( m_fileList.end() - m_deleteBack, m_fileList.end() );
825 m_filesDeleted =
true;
828 if( m_fileList.size() == 0 )
831 "HCIobservation<realT>::readFiles",
832 "The target fileList has 0 length, there are no files to be read after deletions." );
835 if( m_qualityFile !=
"" )
837 std::cerr <<
"Thresholding target images...";
838 size_t origsize = m_fileList.size();
840 if( threshold( m_fileList, m_qualityFile, m_qualityThreshold ) < 0 )
843 if( m_fileList.size() == 0 )
846 "HCIobservation<realT>::readFiles",
847 "The fileList has 0 length, there are no files to be read after thresholding." );
850 std::cerr <<
"Done. Selected " << m_fileList.size() <<
" out of " << origsize <<
"\n";
852 if( m_thresholdOnly )
854 std::cout <<
"#Files which passed thresholding:\n";
855 for(
size_t i = 0; i < m_fileList.size(); ++i )
857 std::cout << m_fileList[i] <<
"\n";
864 if( m_weightFile !=
"" )
866 if( readWeights() < 0 )
874 if( m_MJDKeyword !=
"" )
875 head.
append( m_MJDKeyword );
877 for(
size_t i = 0; i < m_keywords.size(); ++i )
879 if( head.
count( m_keywords[i] ) == 0 )
880 head.
append( m_keywords[i] );
884 m_heads.resize( m_fileList.size(), head );
888 Eigen::Array<realT, Eigen::Dynamic, Eigen::Dynamic> im;
897 m_imSize = im.rows();
898 if( m_imSize > im.cols() )
900 m_imSize = im.cols();
906 if( m_imSize > im.rows() )
908 m_imSize = im.rows();
910 if( m_imSize > im.cols() )
912 m_imSize = im.cols();
918 f.
setReadSize( floor( 0.5 * ( im.rows() - 1 ) - 0.5 * ( m_imSize - 1 ) + 0.1 ),
919 floor( 0.5 * ( im.cols() - 1.0 ) - 0.5 * ( m_imSize - 1.0 ) + 0.1 ),
923 im.resize( m_imSize, m_imSize );
925 std::cerr <<
"image size is " << m_imSize <<
"\n";
931 m_tgtIms.resize( im.rows(), im.cols(), m_fileList.size() );
935 f.
read( m_tgtIms.data(), m_heads, m_fileList );
940 if( m_MJDKeyword !=
"" )
942 m_imageMJD.resize( m_heads.size() );
946 for(
size_t i = 0; i < m_imageMJD.size(); ++i )
953 for(
size_t i = 0; i < m_imageMJD.size(); ++i )
955 m_imageMJD[i] = m_heads[i][m_MJDKeyword].template value<realT>() * m_MJDUnits;
962 m_Nims = m_tgtIms.planes();
963 m_Nrows = m_tgtIms.rows();
964 m_Ncols = m_tgtIms.cols();
965 m_Npix = m_tgtIms.rows() * m_tgtIms.cols();
967 std::cerr <<
"loading complete\n";
969 std::cerr <<
"zero-ing NaNs\n";
971 std::cerr <<
"done\n";
974 if( postReadFiles() < 0 )
983 if( !m_skipPreProcess )
986 if( m_preProcess_beforeCoadd )
988 preProcess( m_tgtIms );
993 std::cerr <<
"Coadding target images...\n";
994 coaddImages( m_coaddCombineMethod,
1004 m_Nims = m_tgtIms.planes();
1005 m_Nrows = m_tgtIms.rows();
1006 m_Ncols = m_tgtIms.cols();
1007 m_Npix = m_tgtIms.rows() * m_tgtIms.cols();
1009 std::cerr <<
"number of target images after coadding: " << m_Nims <<
"\n";
1011 if( postCoadd() < 0 )
1013 std::cerr <<
"Post coadd error " << __FILE__ <<
" " << __LINE__ <<
"\n";
1016 std::cerr <<
"Done.\n";
1019 if( m_maskFile !=
"" )
1026 if( !m_preProcess_beforeCoadd )
1028 preProcess( m_tgtIms );
1031 outputPreProcessed();
1039template <
typename _realT>
1045template <
typename _realT>
1052template <
typename _realT>
1057 if( m_Nrows == 0 || m_Ncols == 0 )
1059 mxError(
"HCIobservation", MXE_PARAMNOTSET,
"The target image size must be set before reading RDI files." );
1063 if( m_RDIfileList.size() == 0 )
1065 mxError(
"HCIobservation", MXE_FILENOTFOUND,
"The RDI fileList has 0 length, there are no files to be read." );
1070 if( !m_RDIfilesDeleted )
1072 if( m_RDIdeleteFront > 0 )
1074 m_RDIfileList.erase( m_RDIfileList.begin(), m_RDIfileList.begin() + m_RDIdeleteFront );
1077 if( m_RDIdeleteBack > 0 )
1079 m_RDIfileList.erase( m_RDIfileList.end() - m_RDIdeleteBack, m_RDIfileList.end() );
1081 m_RDIfilesDeleted =
true;
1084 if( m_RDIfileList.size() == 0 )
1086 mxError(
"HCIobservation",
1088 "The RDI fileList has 0 length, there are no files to be read after deletions." );
1092 if( m_RDIqualityFile !=
"" )
1094 std::cerr <<
"Thresholding RDI images...";
1095 size_t origsize = m_RDIfileList.size();
1097 if( threshold( m_RDIfileList, m_RDIqualityFile, m_RDIqualityThreshold ) < 0 )
1100 if( m_RDIfileList.size() == 0 )
1102 mxError(
"HCIobservation",
1104 "The fileList has 0 length, there are no files to be read after thresholding." );
1108 std::cerr <<
"Done. Selected " << m_RDIfileList.size() <<
" out of " << origsize <<
"\n";
1115 if( m_MJDKeyword !=
"" )
1116 head.
append( m_MJDKeyword );
1118 for(
size_t i = 0; i < m_RDIkeywords.size(); ++i )
1120 head.
append( m_RDIkeywords[i] );
1124 m_RDIheads.resize( m_RDIfileList.size(), head );
1127 Eigen::Array<realT, Eigen::Dynamic, Eigen::Dynamic> im;
1133 if( im.rows() < m_imSize || im.cols() < m_imSize )
1135 mxError(
"HCIobservation", MXE_SIZEERR,
"The reference images are too small, do not match the target images." );
1141 f.
setReadSize( floor( 0.5 * ( im.rows() - 1 ) - 0.5 * ( m_imSize - 1 ) + 0.1 ),
1142 floor( 0.5 * ( im.cols() - 1.0 ) - 0.5 * ( m_imSize - 1.0 ) + 0.1 ),
1146 m_refIms.resize( m_imSize, m_imSize, m_RDIfileList.size() );
1150 f.
read( m_refIms.data(), m_RDIheads, m_RDIfileList );
1154 if( m_MJDKeyword !=
"" )
1156 m_RDIimageMJD.resize( m_RDIheads.size() );
1158 if( m_MJDisISO8601 )
1160 for(
size_t i = 0; i < m_RDIimageMJD.size(); ++i )
1167 for(
size_t i = 0; i < m_RDIimageMJD.size(); ++i )
1169 m_RDIimageMJD[i] = m_RDIheads[i][m_MJDKeyword].template value<realT>() * m_MJDUnits;
1176 std::cerr <<
"loading complete\n";
1178 std::cerr <<
"zero-ing NaNs\n";
1180 std::cerr <<
"Done.\n";
1183 if( postRDIReadFiles() < 0 )
1187 if( !m_skipPreProcess )
1190 if( m_preProcess_beforeCoadd )
1192 preProcess( m_refIms );
1197 std::cerr <<
"Coadding reference images...\n";
1198 coaddImages( m_coaddCombineMethod,
1206 std::cerr <<
"number of reference images after coadding: " << m_refIms.planes() <<
"\n";
1208 if( postRDICoadd() < 0 )
1210 std::cerr <<
"Post coadd error " << __FILE__ <<
" " << __LINE__ <<
"\n";
1213 std::cerr <<
"Done.\n";
1217 if( !m_preProcess_beforeCoadd )
1219 preProcess( m_refIms );
1225 m_RDIfilesRead =
true;
1230template <
typename _realT>
1236template <
typename _realT>
1242template <
typename _realT>
1244 const std::string &qualityFile,
1245 realT qualityThreshold )
1247 if( qualityFile ==
"" )
1249 mxError(
"HCIobservation::threshold", MXE_PARAMNOTSET,
"qualityFile not set" );
1253 int origsize = fileList.size();
1255 std::vector<std::string> qfileNames;
1256 std::vector<realT> imQ;
1261 std::map<std::string, realT> quality;
1262 for(
size_t i = 0; i < qfileNames.size(); ++i )
1267 for(
size_t i = 0; i < fileList.size(); ++i )
1275 q = qualityThreshold - 1;
1278 if( q < qualityThreshold )
1280 fileList.erase( fileList.begin() + i );
1288template <
typename _realT>
1292 std::vector<std::string> &coaddKeywords,
1293 std::vector<double> &imageMJD,
1294 std::vector<fits::fitsHeader> &heads,
1297 std::cerr <<
"***************************************************************\n";
1298 std::cerr <<
" *** WARNING *** \n";
1299 std::cerr <<
" coadding is poorly tested. proceed with caution. \n";
1300 std::cerr <<
"***************************************************************\n";
1303 if( coaddMaxImno <= 0 && coaddMaxTime <= 0 )
1310 int Nims = ims.planes();
1311 int Nrows = ims.rows();
1312 int Ncols = ims.cols();
1316 std::vector<eigenImageT> coadds;
1317 std::vector<std::vector<std::string>> coaddFileNames;
1320 std::vector<double> avgMJD;
1321 std::vector<double> startMJD;
1322 std::vector<double> endMJD;
1324 std::vector<std::vector<double>> avgVals;
1325 std::vector<std::vector<double>> startVals;
1326 std::vector<std::vector<double>> endVals;
1341 std::vector<std::string> imsCoadded;
1349 std::vector<double> initVals;
1350 initVals.resize( coaddKeywords.size() );
1352 std::vector<double> startVal;
1353 startVal.resize( coaddKeywords.size() );
1355 std::vector<double> endVal;
1356 endVal.resize( coaddKeywords.size() );
1361 initMJD = imageMJD[im0];
1362 startMJD.push_back( initMJD );
1363 endMJD.push_back( initMJD );
1365 for(
size_t i = 0; i < coaddKeywords.size(); ++i )
1367 initVals[i] = heads[im0][coaddKeywords[i]].value<
double>();
1368 startVal[i] = initVals[i];
1369 endVal[i] = initVals[i];
1373 bool increment =
true;
1374 while( increment ==
true )
1385 if( imF - im0 > coaddMaxImno && coaddMaxImno > 0 )
1393 if( ( imageMJD[imF] - imageMJD[im0] ) * 86400.0 > coaddMaxTime && coaddMaxTime > 0 )
1405 for(
int imno = im0 + 1; imno < imF; ++imno )
1407 initMJD += imageMJD[imno];
1408 endMJD.back() = imageMJD[imno];
1410 for(
size_t i = 0; i < coaddKeywords.size(); ++i )
1413 heads[imno][coaddKeywords[i]].value<
double>();
1414 initVals[i] += endVal[i];
1419 initMJD /= ( imF - im0 );
1420 for(
size_t i = 0; i < coaddKeywords.size(); ++i )
1422 initVals[i] /= ( imF - im0 );
1426 imsToCoadd.resize( Nrows, Ncols, imF - im0 );
1427 imsCoadded.resize( imF - im0 );
1428 for(
int i = 0; i < ( imF - im0 ); ++i )
1430 imsToCoadd.
image( i ) = ims.
image( im0 + i );
1437 imsToCoadd.
median( coadd );
1441 imsToCoadd.
mean( coadd );
1443 coadds.push_back( coadd );
1444 coaddFileNames.push_back( imsCoadded );
1447 avgMJD.push_back( initMJD );
1448 avgVals.push_back( initVals );
1449 startVals.push_back( startVal );
1450 endVals.push_back( endVal );
1457 ims.resize( Nrows, Ncols, coadds.size() );
1458 Nims = coadds.size();
1460 for(
int i = 0; i < Nims; ++i )
1462 ims.
image( i ) = coadds[i];
1466 imageMJD.erase( imageMJD.begin() + Nims, imageMJD.end() );
1467 heads.erase( heads.begin() + Nims, heads.end() );
1469 for(
int i = 0; i < Nims; ++i )
1471 imageMJD[i] = avgMJD[i];
1475 heads[i].append(
"DELTA " + m_MJDKeyword,
1476 ( endMJD[i] - startMJD[i] ) * 86400,
1477 "change in " + m_MJDKeyword +
" in seconds." );
1479 for(
size_t j = 0; j < coaddKeywords.size(); ++j )
1481 heads[i][coaddKeywords[j]].value( avgVals[i][j] );
1482 heads[i][
"START " + coaddKeywords[j]].value( startVals[i][j] );
1483 heads[i][
"END " + coaddKeywords[j]].value( endVals[i][j] );
1484 heads[i][
"DELTA " + coaddKeywords[j]].value( endVals[i][j] - startVals[i][j] );
1487 heads[i].append(
"IMAGES COADDED", coaddFileNames[i].size(),
"number of images included in this coadd" );
1488 for(
size_t j = 0; j < coaddFileNames[i].size(); ++j )
1492 heads[i].append(
"HISTORY", fits::fitsHistoryType(),
"coadded file: " + coaddFileNames[i][j] );
1500template <
typename _realT>
1504 if( m_maskFile !=
"" )
1506 std::cerr <<
"creating mask cube\n";
1508 ff.
read( m_mask, m_maskFile );
1511 if( m_mask.rows() > m_imSize || m_mask.cols() > m_imSize )
1513 eigenImageT tmask = m_mask.block( (
int)( 0.5 * ( m_mask.rows() - 1 ) - 0.5 * ( m_imSize - 1 ) ),
1514 (
int)( 0.5 * ( m_mask.rows() - 1 ) - 0.5 * ( m_imSize - 1 ) ),
1524template <
typename realT>
1527 if( m_mask.rows() != m_Nrows || m_mask.cols() != m_Ncols )
1530 std::string message =
"Mask is not the same size as images.\n";
1531 message +=
" Mask: rows=" + std::to_string(m_mask.rows()) +
"\n";
1532 message +=
" cols=" + std::to_string(m_mask.cols()) +
"\n";
1533 message +=
" Images: rows=" + std::to_string(m_Nrows) +
"\n";
1534 message +=
" cols=" + std::to_string(m_Ncols) +
"\n";
1537 mxThrowException(
err::invalidconfig,
"HCIobservation<realT>::makeMaskCube", message );
1540 m_maskCube.resize( m_Nrows, m_Ncols, m_Nims );
1542 for(
int i = 0; i < m_Nims; ++i )
1544 m_maskCube.image( i ) = m_mask;
1548template <
typename _realT>
1554 if( m_maskFile !=
"" && m_preProcess_mask )
1556 std::cerr <<
"Masking . . .\n";
1557#pragma omp parallel for
1558 for(
int i = 0; i < ims.planes(); ++i )
1560 ims.
image( i ) *= m_mask;
1562 std::cerr <<
"done\n";
1565 if( m_preProcess_subradprof )
1567 std::cerr <<
"subtracting radial profile . . .\n";
1574 for(
int i = 0; i < ims.planes(); ++i )
1576 Eigen::Map<Eigen::Array<realT, Eigen::Dynamic, Eigen::Dynamic>> imRef( ims.
image( i ) );
1582 std::cerr <<
"done\n";
1584 if( m_maskFile !=
"" && m_preProcess_mask )
1586 std::cerr <<
"Masking . . .\n";
1587#pragma omp parallel for
1588 for(
int i = 0; i < ims.planes(); ++i )
1590 ims.
image( i ) *= m_mask;
1592 std::cerr <<
"done\n";
1596 if( m_preProcess_medianUSM_fwhm > 0 )
1598 std::cerr <<
"applying median USM . . .\n";
1601 medmask.resize( ims.rows(), ims.cols() );
1602 medmask.setConstant( 1 );
1604 for(
int z = 0; z < (int)( 0.5 * m_preProcess_medianUSM_fwhm ); ++z )
1606 medmask.row( z ) = 0;
1607 medmask.row( medmask.rows() - 1 - z ) = 0;
1608 medmask.col( z ) = 0;
1609 medmask.col( medmask.cols() - 1 - z ) = 0;
1615 fim.resize( ims.rows(), ims.cols() );
1616 im.resize( ims.rows(), ims.cols() );
1619 for(
int i = 0; i < ims.planes(); ++i )
1621 im = ims.
image( i );
1623 ims.
image( i ) = ( im - fim ) * medmask;
1627 if( m_maskFile !=
"" && m_preProcess_mask )
1629 std::cerr <<
"Masking . . .\n";
1630#pragma omp parallel for
1631 for(
int i = 0; i < ims.planes(); ++i )
1633 ims.
image( i ) *= m_mask;
1636 std::cerr <<
"done\n";
1639 if( m_preProcess_gaussUSM_fwhm > 0 )
1641 std::cerr <<
"applying Gauss USM . . .\n";
1644#pragma omp parallel for
1645 for(
int i = 0; i < ims.planes(); ++i )
1648 im = ims.
image( i );
1652 0.5 * ( ims.cols() - 1 ) - m_preProcess_gaussUSM_fwhm * 4 );
1654 ims.
image( i ) = im;
1657 if( m_maskFile !=
"" && m_preProcess_mask )
1659 std::cerr <<
"Masking . . .\n";
1660#pragma omp parallel for
1661 for(
int i = 0; i < ims.planes(); ++i )
1663 ims.
image( i ) *= m_mask;
1667 std::cerr <<
"done\n";
1670 if( m_preProcess_azUSM_azW && m_preProcess_azUSM_radW )
1674 std::cerr <<
"applying azimuthal USM . . .\n";
1676#pragma omp parallel for
1677 for(
int i = 0; i < ims.planes(); ++i )
1680 im = ims.
image( i );
1681 medianFilterImage( fim,
1684 m_preProcess_azUSM_azW,
1685 m_preProcess_azUSM_maxAz ) );
1687 ims.
image( i ) = im;
1691 if( m_maskFile !=
"" && m_preProcess_mask )
1693 std::cerr <<
"Masking . . .\n";
1694#pragma omp parallel for
1695 for(
int i = 0; i < ims.planes(); ++i )
1697 ims.
image( i ) *= m_mask;
1702 std::cerr <<
"done (" << t_azusm_end - t_azusm_begin <<
" sec) \n";
1705 preProcess_meanSub( ims );
1707 preProcess_pixelTSNorm( ims );
1713template <
typename _realT>
1723 std::string msg =
"Mean subtraction by " + HCI::meanSubMethodStr( m_preProcess_meanSubMethod );
1724 msg +=
" can't be done in pre-processing. Only meanImage or medianImage can be used in pre.";
1739#pragma omp parallel for
1740 for(
int n = 0; n < ims.planes(); ++n )
1742 ims.
image( n ) -= mean;
1744 if( m_maskFile !=
"" && m_preProcess_mask )
1746 ims.
image( n ) *= m_mask;
1751template <
typename _realT>
1761 mxThrowException(
err::invalidconfig,
"KlipReduction::preProcess_pixelTSNorm",
"pixelTSNormMethod is unknown" );
1767 "KlipReduction::preProcess_pixelTSNorm",
1768 "pixelTSNormMethod is rmsSigmaClipped, which is not implemented" );
1771 std::cerr <<
"normalizing pixels\n";
1775 std::vector<realT> pixs( ims.planes() );
1778 for(
int cc = 0; cc < ims.cols(); ++cc )
1780 for(
int rr = 0; rr < ims.rows(); ++rr )
1782 if( m_maskFile !=
"" && m_preProcess_mask )
1784 if( m_mask( rr, cc ) == 0 )
1791 for(
int pp = 0; pp < ims.planes(); ++pp )
1793 pixs[pp] = ims.
image( pp )( rr, cc );
1801 for(
int pp = 0; pp < ims.planes(); ++pp )
1803 ims.
image( pp )( rr, cc ) /= sd;
1810template <
typename _realT>
1816 if( m_weightFile ==
"" )
1818 mxError(
"HCIobservation::readWeights", MXE_PARAMNOTSET,
"m_weightFile not set" );
1823 std::vector<std::string> wfileNames;
1824 std::vector<realT> imW;
1829 if( imW.size() < m_fileList.size() )
1831 mxError(
"HCIobservation::readWeights", MXE_SIZEERR,
"not enough weights specified" );
1835 std::map<std::string, realT> weights;
1836 for(
size_t i = 0; i < wfileNames.size(); ++i )
1839 m_comboWeights.resize( m_fileList.size() );
1842 realT weightSum = 0;
1843 for(
size_t i = 0; i < m_fileList.size(); ++i )
1851 mxError(
"HCIobservation::readWeights", MXE_NOTFOUND,
"Weight for a file in m_fileList not found." );
1854 m_comboWeights[i] = wi;
1859 for(
size_t i = 0; i < m_comboWeights.size(); ++i )
1861 m_comboWeights[i] /= weightSum;
1867template <
typename _realT>
1890 if( m_sigmaThreshold > 0 )
1899 else if( m_combineMethod == HCI::debug )
1901 method = HCI::debug;
1907 m_finim.resize( m_psfsub[0].rows(), m_psfsub[0].cols(), m_psfsub.size() );
1910 for(
size_t n = 0; n < m_psfsub.size(); ++n )
1914 m_psfsub[n].median( tfinim );
1915 m_finim.image( n ) = tfinim;
1919 if( m_comboWeights.size() == (
size_t)m_Nims )
1921 if( m_maskFile !=
"" )
1923 m_psfsub[n].mean( tfinim, m_comboWeights, m_maskCube, m_minGoodFract );
1927 m_psfsub[n].mean( tfinim, m_comboWeights );
1932 if( m_maskFile !=
"" )
1934 m_psfsub[n].mean( tfinim, m_maskCube, m_minGoodFract );
1938 m_psfsub[n].mean( tfinim );
1941 m_finim.image( n ) = tfinim;
1945 if( m_comboWeights.size() == (
size_t)m_Nims )
1947 if( m_maskFile !=
"" )
1949 m_psfsub[n].sigmaMean( tfinim, m_comboWeights, m_maskCube, m_sigmaThreshold, m_minGoodFract );
1953 m_psfsub[n].sigmaMean( tfinim, m_comboWeights, m_sigmaThreshold );
1958 if( m_maskFile !=
"" )
1960 m_psfsub[n].sigmaMean( tfinim, m_maskCube, m_sigmaThreshold, m_minGoodFract );
1964 m_psfsub[n].sigmaMean( tfinim, m_sigmaThreshold );
1967 m_finim.image( n ) = tfinim;
1969 else if( method == HCI::debug )
1971 m_finim.image( n ) = m_psfsub[n].image( 0 );
1978template <
typename _realT>
1981 if( m_preProcess_outputPrefix ==
"" )
1995 for(
int i = 0; i < m_Nims; ++i )
1997 snprintf( nstr,
sizeof( nstr ),
"%06d", i );
1998 fname = m_preProcess_outputPrefix + nstr +
".fits";
2001 stdFitsHeader( fh );
2002 ff.
write( fname, m_tgtIms.image( i ).data(), m_Ncols, m_Nrows, 1, fh );
2006template <
typename _realT>
2009 head.
append(
"", fits::fitsCommentType(),
"----------------------------------------" );
2010 head.
append(
"", fits::fitsCommentType(),
"mx::HCIobservation parameters:" );
2011 head.
append(
"", fits::fitsCommentType(),
"----------------------------------------" );
2013 head.
append<
int>(
"FDELFRNT", m_deleteFront,
"images deleted from front of file list" );
2014 head.
append<
int>(
"FDELBACK", m_deleteBack,
"images deleted from back of file list" );
2017 head.
append<
realT>(
"QTHRESH", m_qualityThreshold,
"quality threshold" );
2018 head.
append<
int>(
"NUMIMS", m_Nims,
"number of images processed" );
2020 head.
append<
int>(
"IMSIZE", m_imSize,
"image size after reading" );
2025 head.
append<
int>(
"COADIMNO", m_coaddMaxImno,
"max number of images in each coadd" );
2026 head.
append<
realT>(
"COADTIME", m_coaddMaxTime,
"max time in each coadd" );
2030 head.
append<
int>(
"COADIMNO", 0,
"max number of images in each coadd" );
2031 head.
append<
realT>(
"COADTIME", 0,
"max time in each coadd" );
2034 head.
append(
"MASKFILE", m_maskFile,
"mask file" );
2036 head.
append<
int>(
"PREPROC BEFORE", m_preProcess_beforeCoadd,
"pre-process before coadd flag" );
2037 head.
append<
int>(
"PREPROC MASK", m_preProcess_mask,
"pre-process mask flag" );
2038 head.
append<
int>(
"PREPROC SUBRADPROF", m_preProcess_subradprof,
"pre-process subtract radial profile flag" );
2040 m_preProcess_azUSM_azW,
2041 "pre-process azimuthal USM azimuthal width [pixels]" );
2043 m_preProcess_azUSM_maxAz,
2044 "pre-process azimuthal USM maximum azimuthal width [degrees]" );
2046 m_preProcess_azUSM_radW,
2047 "pre-process azimuthal USM radial width [pixels]" );
2048 head.
append<
realT>(
"PREPROC MEDIANUSM FWHM", m_preProcess_medianUSM_fwhm,
"pre-process median USM fwhm [pixels]" );
2049 head.
append<
realT>(
"PREPROC GAUSSUSM FWHM", m_preProcess_gaussUSM_fwhm,
"pre-process Gaussian USM fwhm [pixels]" );
2050 head.
append<std::string>(
"PREPROC MEANSUB METHOD",
2051 HCI::meanSubMethodStr( m_preProcess_meanSubMethod ),
2052 "pre-process mean subtraction method" );
2053 head.
append<std::string>(
"PREPROC PIXELTSNORM METHOD",
2054 HCI::pixelTSNormMethodStr( m_preProcess_pixelTSNormMethod ),
2055 "pre-process pixel time-series norm method" );
2058template <
typename _realT>
2061 std::string fname = m_finimName;
2063 if( m_outputDir !=
"" )
2065 mkdir( m_outputDir.c_str(), S_IRWXU | S_IRWXG | S_IROTH | S_IXOTH );
2067 fname = m_outputDir +
"/" + fname;
2070 if( !m_exactFinimName )
2078 stdFitsHeader( head );
2083 if( m_weightFile !=
"" )
2084 head.
append(
"WEIGHT FILE", m_weightFile,
"file containing weights for combination" );
2087 head.
append<
realT>(
"SIGMA THRESHOLD", m_sigmaThreshold,
"threshold for sigma clipping" );
2089 head.
append<
realT>(
"MIN FOOD FRACTION", m_minGoodFract,
"minimum good fraction for inclusion" );
2099 f.
write( fname, m_finim, head );
2101 std::cerr <<
"Final image written to: " << fname <<
"\n";
2104template <
typename _realT>
2113 stdFitsHeader( head );
2126 if( m_comboWeights.size() > 0 )
2128 fname = m_PSFSubPrefix +
"weights.dat";
2129 if( m_outputDir !=
"" )
2131 mkdir( m_outputDir.c_str(), S_IRWXU | S_IRWXG | S_IROTH | S_IXOTH );
2132 fname = m_outputDir +
"/" + fname;
2135 std::cerr <<
"writing comboWeights: " << fname <<
"\n";
2139 for(
size_t n = 0; n < m_psfsub.size(); ++n )
2141 for(
int p = 0; p < m_psfsub[n].planes(); ++p )
2143 snprintf( num, 256,
"_%03zu_%05d.fits", n, p );
2144 fname = m_PSFSubPrefix + num;
2146 if( m_outputDir !=
"" )
2148 mkdir( m_outputDir.c_str(), S_IRWXU | S_IRWXG | S_IROTH | S_IXOTH );
2149 fname = m_outputDir +
"/" + fname;
2154 h.append( m_heads[p] );
2156 f.
write( fname, m_psfsub[n].image( p ).data(), m_psfsub[n].rows(), m_psfsub[n].cols(), 1, h );
2158 if( m_comboWeights.size() > 0 && n == 0 )
2160 wout << fname <<
" " << m_comboWeights[p] <<
"\n";
2165 if( m_comboWeights.size() > 0 )
2171template <
typename _realT>
2173 const std::string &prefix,
2174 const std::string &ext,
2175 size_t nReductions )
2178 m_psfsub.resize( nReductions );
2186 ff.
read( im, fh, flist[0] );
2188 if( fh.
count(
"FDELFRNT" ) == 0 )
2190 mxError(
"KLIPReduction", MXE_PARAMNOTSET,
"FDELFRNT not found in FITS header." );
2193 m_deleteFront = fh[
"FDELFRNT"].value<
int>();
2194 std::cerr <<
"deleteFront: " << m_deleteFront <<
"\n";
2196 if( fh.
count(
"FDELBACK" ) == 0 )
2198 mxError(
"KLIPReduction", MXE_PARAMNOTSET,
"FDELBACK not found in FITS header." );
2201 m_deleteBack = fh[
"FDELBACK"].value<
int>();
2202 std::cerr <<
"deleteBack: " << m_deleteBack <<
"\n";
2204 if( fh.
count(
"QFILE" ) == 0 )
2206 mxError(
"KLIPReduction", MXE_PARAMNOTSET,
"QFILE not found in FITS header." );
2209 m_qualityFile = fh[
"QFILE"].String();
2210 std::cerr <<
"qualityFile: " << m_qualityFile <<
"\n";
2212 if( fh.
count(
"QTHRESH" ) == 0 )
2214 mxError(
"KLIPReduction", MXE_PARAMNOTSET,
"QTHRESH not found in FITS header." );
2217 m_qualityThreshold = fh[
"QTHRESH"].value<
realT>();
2218 std::cerr <<
"qualityThreshold: " << m_qualityThreshold <<
"\n";
2220 if( fh.
count(
"COADMTHD" ) == 0 )
2222 mxError(
"KLIPReduction", MXE_PARAMNOTSET,
"COADMTHD not found in FITS header." );
2226 std::cerr <<
"coaddCombineMethod: " << m_coaddCombineMethod <<
"\n";
2228 if( fh.
count(
"COADIMNO" ) != 0 )
2230 m_coaddMaxImno = fh[
"COADIMNO"].value<
int>();
2231 std::cerr <<
"coaddMaxImno: " << m_coaddMaxImno <<
"\n";
2234 if( fh.
count(
"COADTIME" ) != 0 )
2236 m_coaddMaxImno = fh[
"COADTIME"].value<
realT>();
2237 std::cerr <<
"coaddMaxtime: " << m_coaddMaxTime <<
"\n";
2240 if( m_maskFile ==
"" )
2242 if( fh.
count(
"MASKFILE" ) == 0 )
2244 mxError(
"KLIPReduction",
2246 "MASKFILE not found in FITS header and not set in configuration." );
2249 m_maskFile = fh[
"MASKFILE"].String();
2251 std::cerr <<
"maskFile: " << m_maskFile <<
"\n";
2253 if( fh.
count(
"PPBEFORE" ) == 0 )
2255 mxError(
"KLIPReduction", MXE_PARAMNOTSET,
"PPBEFORE not found in FITS header." );
2258 m_preProcess_beforeCoadd = fh[
"PPBEFORE"].value<
int>();
2259 std::cerr <<
"preProcess_beforeCoadd: " << m_preProcess_beforeCoadd <<
"\n";
2261 if( fh.
count(
"PPMASK" ) == 0 )
2263 mxError(
"KLIPReduction", MXE_PARAMNOTSET,
"PPMASK not found in FITS header." );
2266 m_preProcess_mask = fh[
"PPMASK"].value<
int>();
2267 std::cerr <<
"preProcess_mask: " << m_preProcess_mask <<
"\n";
2269 if( fh.
count(
"PPSUBRAD" ) == 0 )
2271 mxError(
"KLIPReduction", MXE_PARAMNOTSET,
"PPSUBRAD not found in FITS header." );
2274 m_preProcess_subradprof = fh[
"PPSUBRAD"].value<
int>();
2275 std::cerr <<
"preProcess_subradprof: " << m_preProcess_subradprof <<
"\n";
2277 if( fh.
count(
"PPAUSMAW" ) == 0 )
2279 mxError(
"KLIPReduction", MXE_PARAMNOTSET,
"PPAUSMAW not found in FITS header." );
2282 m_preProcess_azUSM_azW = fh[
"PPAUSMAW"].value<
realT>();
2283 std::cerr <<
"preProcess_azUSM_azW: " << m_preProcess_azUSM_azW <<
"\n";
2285 if( fh.
count(
"PPAUSMRW" ) == 0 )
2287 mxError(
"KLIPReduction", MXE_PARAMNOTSET,
"PPAUSMRW not found in FITS header." );
2290 m_preProcess_azUSM_radW = fh[
"PPAUSMRW"].value<
realT>();
2291 std::cerr <<
"preProcess_azUSM_radW: " << m_preProcess_azUSM_radW <<
"\n";
2293 if( fh.
count(
"PPGUSMFW" ) == 0 )
2295 mxError(
"KLIPReduction", MXE_PARAMNOTSET,
"PPGUSMFW not found in FITS header." );
2298 m_preProcess_gaussUSM_fwhm = fh[
"PPGUSMFW"].value<
realT>();
2299 std::cerr <<
"preProcess_gaussUSM_fwhm: " << m_preProcess_gaussUSM_fwhm <<
"\n";
2303 if( m_MJDKeyword !=
"" )
2304 head.
append( m_MJDKeyword );
2306 for(
size_t i = 0; i < m_keywords.size(); ++i )
2308 head.
append( m_keywords[i] );
2311 for(
size_t n = 0; n < nReductions; ++n )
2314 int nwr = snprintf( nstr,
sizeof( nstr ),
"%03zu", n );
2315 if( nwr < 0 || n >=
sizeof( nstr ) )
2317 std::cerr <<
"possibly bad formatting in filename\n";
2320 std::string nprefix = prefix +
"_" + nstr +
"_";
2321 load_fileList( dir, nprefix, ext );
2323 if( m_fileList.size() == 0 )
2325 mxError(
"HCIobservation",
2327 "The m_fileList has 0 length, there are no files to be read." );
2331 Eigen::Array<realT, Eigen::Dynamic, Eigen::Dynamic> im;
2341 m_imSize = im.rows();
2342 if( m_imSize > im.cols() )
2344 m_imSize = im.cols();
2350 if( m_imSize > im.rows() )
2352 m_imSize = im.rows();
2354 if( m_imSize > im.cols() )
2356 m_imSize = im.cols();
2361 f.
setReadSize( floor( 0.5 * ( im.rows() - 1 ) - 0.5 * ( m_imSize - 1 ) + 0.1 ),
2362 floor( 0.5 * ( im.cols() - 1.0 ) - 0.5 * ( m_imSize - 1.0 ) + 0.1 ),
2365 im.resize( m_imSize, m_imSize );
2367 std::cerr <<
"image size is " << m_imSize <<
"\n";
2371 if( m_fileList.size() != (
size_t)m_Nims )
2373 mxError(
"HCIobservation", MXE_INVALIDARG,
"Different number of images in reductions." );
2376 if( m_Nrows != im.rows() )
2378 mxError(
"HCIobservation", MXE_INVALIDARG,
"Different number of rows in reductions." );
2381 if( m_Ncols != im.cols() )
2383 mxError(
"HCIobservation", MXE_INVALIDARG,
"Different number of cols in reductions." );
2389 std::cerr <<
"found " << nReductions <<
" sets of " << m_fileList.size() <<
" " << im.rows() <<
" x "
2390 << im.cols() <<
" files\n";
2392 m_Nims = m_fileList.size();
2393 m_Nrows = im.rows();
2394 m_Ncols = im.cols();
2395 m_Npix = im.rows() * im.cols();
2397 m_psfsub[n].resize( m_Nrows, m_Ncols, m_Nims );
2400 m_heads.resize( m_fileList.size(), head );
2404 f.
read( m_psfsub[n].data(), m_heads, m_fileList );
2408 if( m_MJDKeyword !=
"" )
2410 m_imageMJD.resize( m_heads.size() );
2412 if( m_MJDisISO8601 )
2414 for(
size_t i = 0; i < m_imageMJD.size(); ++i )
2421 for(
size_t i = 0; i < m_imageMJD.size(); ++i )
2423 m_imageMJD[i] = m_heads[i][m_MJDKeyword].template value<realT>() * m_MJDUnits;
2430 for(
size_t n = 0; n < m_psfsub.size(); ++n )
2436 if( m_weightFile !=
"" )
2438 std::vector<std::string> fn;
2441 std::cerr <<
"read: " << m_weightFile <<
" (" << m_comboWeights.size() <<
")\n";
2445 if( postReadFiles() < 0 )
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.
void setReadSize()
Set to read all the pixels in the 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.
void mean(eigenT &mim)
Calculate the mean image of the cube.
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.
An image cube with an Eigen API.
Tools for using the eigen library for image processing.
int readColumns(const std::string &fname, arrTs &...arrays)
Read in columns from a text file.
Eigen::Array< scalarT, -1, -1 > eigenImage
Definition of the eigenImage type, which is an alias for Eigen::Array.
std::string getSequentialFilename(const std::string &basename, const std::string &extension="", const int startat=0, int ndigit=4)
Get the next file in a numbered sequence.
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,...
int createDirectories(const std::string &path)
Create a directory or directories.
std::string parentPath(const std::string &fname)
Get the parent path from a filename.
std::string pathFilename(const std::string &fname)
Get the base filename.
void fitsHeaderGitStatus(fitsHeader &head, const std::string &repoName, const char *sha1, int modified)
Write the status of a Git repository to HISTORY in a FITS header.
combineMethods
Possible combination methods.
meanSubMethod
Mean subtraction methods.
@ sigmaMeanCombine
Combine with the sigma clipped mean.
@ noCombine
Do not combine the images.
@ meanCombine
Combine with the mean.
@ medianCombine
Combine with the median.
@ none
No mean subtraction.
@ medianImage
The median image of the data is subtracted from each image.
@ meanImage
The mean image of the data is subtracted from each image.
@ unknown
unknown value, an error
int medianSmooth(imageTout &imOut, int &xMax, int &yMax, typename imageTout::Scalar &pMax, const imageTin &imIn, int medianFullWidth)
void filterImage(imageOutT &fim, imageInT im, kernelT kernel, int maxr=0)
Filter an image with a mean kernel.
void zeroNaNs(imageT &im, valueT val)
Zero any NaNs in an image.
void zeroNaNCube(cubeT &imc, maskCubeT *mask)
Zero any NaNs in an image cube.
void radprofim(radprofT &radprofIm, eigenImT1 &im, const eigenImT2 &rad, const eigenImT3 *mask, bool subtract, bool mean=false)
Form a radial profile image, and optionally subtract it from the input.
std::string ISO8601DateTimeStrMJD(const double &timeIn, int timeZone=0)
Get a date-time string in ISO 8601 format for an MJD.
typeT get_curr_time()
Get the current system time in seconds.
double ISO8601date2mjd(const std::string &fdate)
Parse an ISO8601 date of the form "YYYY-MM-DDTHH:MM:SS.S" and return the modified Julian date (MJD)
valueT vectorVariance(const valueT *vec, size_t sz, valueT mean)
Calculate the variance of a vector relative to a supplied mean value.
Image filters (smoothing, radial profiles, etc.)
Declares and defines functions to work with image masks.
Header for the image processing utilities.
int combineMethodFmStr(const std::string &method)
Get the combineMethod from the corresponding string name.
std::string combineMethodStr(int method)
Get the string name of the combineMethod integer.
@ rmsSigmaClipped
the sigma clipped rms of the pixel time series
@ none
no pixel time series norm
@ unknown
unknown value, an error
@ rms
the rms of the pixel time series
The basic high contrast imaging data type.
int m_imSize
Set the image size. Images are cut down to this size after reading.
int m_RDIdeleteBack
Specify how many files from m_RDIfileList to delete from the back of the list.
bool m_preProcess_only
If true, then we stop after pre-processing.
_realT realT
The arithmetic type used for calculations. Does not have to match the type in images on disk.
std::vector< std::string > m_coaddKeywords
The values of these keywords will be averaged and replaced.
bool m_thresholdOnly
Just prints the names and qualities of the files which pass threshold, and stop.
void stdFitsHeader(fits::fitsHeader &head)
Fill in the HCIobservation standard FITS header.
std::string m_weightFile
Specifies a file containing the image weights, for combining with weighted mean.
std::vector< std::string > m_RDIfileList
The list of files to read in for the RDI basis.
std::vector< std::string > m_RDIkeywords
Vector of FITS header keywords to read from the files in m_fileList.
bool m_preProcess_mask
If true, the mask is applied during each pre-processing step.
std::string m_RDIqualityFile
File containing 2 space-delimited columns of fileMame qualityValue pairs for the reference images.
eigenCube< realT > m_tgtIms
The target image cube.
bool m_RDIfilesRead
Whether or not the reference files have been read.
std::vector< double > m_RDIimageMJD
Vector of reference image times, in MJD.
HCIobservation()
Default c'tor.
realT m_MJDUnits
If the date is not ISO 8601, this specifies the conversion to Julian Days (i.e. seconds to days)
int m_preProcess_medianUSM_fwhm
Kernel FWHM for symmetric box median unsharp mask (USM)
bool m_exactFinimName
Use m_finimName exactly as specified, without appending a number or an extension.
int readRDIFiles()
Read the list of reference files, cut to size, and preprocess.
void writeFinim(fits::fitsHeader *addHead=0)
Write the final combined image to disk.
std::vector< fits::fitsHeader > m_RDIheads
Vector of FITS headers, one per reference file, populated with the values for the keywords.
std::string m_preProcess_outputPrefix
Set path and file prefix to output the pre-processed images.
int m_Nrows
Number of rows of the images in m_tgtIms.
void load_RDIfileList(const std::string &dir, const std::string &prefix, const std::string &ext)
Load the RDI basis file list.
int readPSFSub(const std::string &dir, const std::string &prefix, const std::string &ext, size_t nReductions)
Read in already PSF-subtracted files.
realT m_pixelTSSigma
Sigma-clipping parameter for pixel time-series normalization.
int m_deleteBack
Specify how many files from m_fileList to delete from the back of the list.
void preProcess_meanSub(eigenCube< realT > &ims)
Do mean subtraction as part of pre-processing.
HCI::meanSubMethod m_preProcess_meanSubMethod
The mean subtraction method during pre-processing.
std::string m_qualityFile
File containing 2 space-delimited columns of fileVame qualityValue pairs.
eigenImageT m_mask
The mask.
int readWeights()
Read the image weights from m_weightFile.
int m_combineMethod
Determine how to combine the PSF subtracted images.
Eigen::Array< realT, Eigen::Dynamic, Eigen::Dynamic > eigenImageT
The Eigen image array type basted on realT.
realT m_preProcess_azUSM_radW
Radial boxcar width for azimuthal unsharp mask [pixels].
std::string m_auxDataDir
Location for temporary auxilliary output files (e.g. masks)
realT m_sigmaThreshold
The standard deviation threshold used if combineMethod == HCI::sigmaMeanCombine.
std::string m_finimName
The base file name of the output final image.
bool m_skipPreProcess
Don't do any of the pre-processing steps (including coadding).
int m_Ncols
Number of columns of the images in m_tgtIms.
int m_Nims
Number of images in m_tgtIms.
void readMask()
Read the mask file, resizing to imSize if needed.
int m_coaddMaxImno
Maximum number of images to coadd.
bool m_MJDisISO8601
Whether or not the date is in ISO 8601 format.
bool m_filesRead
Whether or not the m_fileList has been read.
std::string m_outputDir
The directory where to write output files.
bool m_preProcess_subradprof
If true, a radial profile is subtracted from each image.
std::vector< std::string > m_fileList
The list of files to read in.
realT m_coaddMaxTime
Maximum elapsed time over which to coadd the images.
eigenCube< realT > m_maskCube
void outputRDIPreProcessed()
Output the pre-processed reference images.
virtual int postCoadd()
Perform post-coadd actions for the target images, for use by derived classes.
HCI::pixelTSNormMethod m_preProcess_pixelTSNormMethod
Specify if each pixel time-series is normalized.
void outputPSFSub(fits::fitsHeader *addHead=0)
Write the PSF subtracted images to disk.
std::vector< eigenCube< realT > > m_psfsub
The PSF subtracted images.
std::string m_PSFSubPrefix
bool m_preProcess_beforeCoadd
controls whether pre-processing takes place before or after coadding
void preProcess_pixelTSNorm(eigenCube< realT > &ims)
Do pixel time-series normalization as part of pre-processing.
void load_fileList(const std::string &dir, const std::string &prefix, const std::string &ext)
Load the file list.
int m_Npix
Pixels per image, that is Nrows*Ncols.
int threshold(std::vector< std::string > &fileList, const std::string &qualityFile, realT qualityThreshold)
Read the image qualities from a qualityFile and apply the threshold to a fileList.
std::string m_MJDKeyword
Name of the keyword to use for the image date.
eigenCube< realT > m_refIms
The optional reference image cube.
std::vector< std::string > m_keywords
Vector of FITS header keywords to read from the files in m_fileList.
bool m_doOutputPSFSub
Controls whether or not the individual PSF subtracted images are written to disk.
realT m_RDIqualityThreshold
Threshold to apply to qualityValues read from qualityFile.
int m_coaddCombineMethod
Determine how to coadd the raw images.
realT m_qualityThreshold
Threshold to apply to qualityValues read from qualityFile.
virtual int postRDIReadFiles()
Perform post-read actions for the RDI images, for use by derived classes.
int m_deleteFront
Specify how many files from m_fileList to delete from the front of the list.
virtual void makeMaskCube()
Populate the mask cube which is used for post-processing.
realT m_preProcess_gaussUSM_fwhm
Kernel FWHM for symmetric Gaussian unsharp mask (USM)
virtual int postReadFiles()
Perform post-read actions for the target images, for use by derived classes.
realT m_preProcess_azUSM_azW
Azimuthal boxcar width for azimuthal unsharp mask [pixels].
void coaddImages(int coaddCombineMethod, int coaddMaxImno, int coaddMaxTime, std::vector< std::string > &coaddKeywords, std::vector< double > &imageMJD, std::vector< fits::fitsHeader > &heads, eigenCube< realT > &ims)
Coadd the images.
std::vector< realT > m_comboWeights
Vector to hold the image weights read from the m_weightFile.
bool m_moveAuxDataDir
Whether or not to move the temp. aux files.
std::vector< fits::fitsHeader > m_heads
Vector of FITS headers, one per file, populated with the values for the keywords.
std::string m_maskFile
Specify a mask file to apply.
int m_doWriteFinim
Set whether the final combined image is written to disk.
virtual int postRDICoadd()
Perform post-coadd actions, for use by derived classes.
bool m_filesDeleted
Whether or not the specified files have been deleted from m_fileList.
void combineFinim()
Combine the images into a single final image.
realT m_preProcess_azUSM_maxAz
Mazimum azimuthal boxcar width for azimuthal unsharp mask [degrees].
void preProcess(eigenCube< realT > &ims)
Do the pre-processing.
int readFiles()
Read the list of files, cut to size, and preprocess.
int m_RDIdeleteFront
Specify how many files from m_RDIfileList to delete from the front of the list.
std::vector< double > m_imageMJD
Vector of target image times, in MJD.
void outputPreProcessed()
Output the pre-processed target images.
eigenCube< realT > m_finim
The final combined images, one for each cube in psfsub.
bool m_RDIfilesDeleted
Whether or not the specified files have been deleted from m_RDIfileList.
Azimuthally variable boxcare kernel.
Symetric Gaussian smoothing kernel.