mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
Loading...
Searching...
No Matches
KLIPreduction.hpp
Go to the documentation of this file.
1/** \file KLIPreduction.hpp
2 * \author Jared R. Males
3 * \brief Declarations and definitions for an implementation of the
4 * Karhunen-Loeve Image Processing (KLIP) algorithm. \ingroup hc_imaging_files
5 * \ingroup image_processing_files
6 *
7 */
8
9#ifndef __KLIPreduction_hpp__
10#define __KLIPreduction_hpp__
11
12#include <map>
13#include <vector>
14
15#include <omp.h>
16
17#include "../ipc/ompLoopWatcher.hpp"
18#include "../math/eigenLapack.hpp"
19#include "../math/geo.hpp"
20#include "../sigproc/gramSchmidt.hpp"
21using namespace mx::sigproc;
22
23#include "ADIobservation.hpp"
24
25namespace mx
26{
27namespace improc
28{
29
30// double t0, t1, t2, t3, t4, t5, t6, t7, t8, t9, t10, t11, tf;
31// double dcut, dscv, dklims, dgemm, dsyevr, dcfs, drot, dcombo, dread;
32
33namespace HCI
34{
35
36
37
38/// Image exclusion methods
39/** \ingroup hc_imaging_enums
40 */
42{
43 excludeNone, ///< Exclude no images
44 excludePixel, ///< Exclude by pixels of rotation
45 excludeAngle, ///< Exclude by angle of roration
46 excludeImno ///< Exclude by number of images
47};
48
49std::string excludeMethodStr( int method );
50
51int excludeMethodFmStr( const std::string &method );
52
53/// Image inclusion methods
54/** \ingroup hc_imaging_enums
55 */
57{
58 includeAll, ///< include all images
59 includeCorr, ///< include images which are most correlated with the target
60 includeTime, ///< include images which are closest in time to the target
61 includeAngle, ///< include images which are closest in angle to the target
62 includeImno ///< include images which are closest in imno to the target
63};
64
65std::string includeMethodStr( int method );
66
67int includeMethodFmStr( const std::string &method );
68} // namespace HCI
69
70
71
72/// An implementation of the Karhunen-Loeve Image Processing (KLIP) algorithm.
73/** KLIP\cite soummer_2012 is a principle components analysis (PCA) based
74 * technique for PSF estimation.
75 *
76 *
77 * \tparam _realT is the floating point type in which to do calculations
78 * \tparam _derotFunctObj the ADIobservation derotator class.
79 * \tparam _evCalcT the real type in which to do eigen-decomposition. Should
80 * generally be double for stable results. \ingroup hc_imaging
81 */
82template <typename _realT, class _derotFunctObj, typename _evCalcT = double>
83struct KLIPreduction : public ADIobservation<_realT, _derotFunctObj>
84{
85 typedef _realT realT;
86
87 typedef Eigen::Array<realT, Eigen::Dynamic, Eigen::Dynamic> eigenImageT;
88
89 typedef _evCalcT evCalcT;
90
91 /// Default c'tor
93
94 /// Construct and load the target file list.
95 /** Populates the \ref m_fileList vector by searching on disk for files
96 * which match "dir/prefix*.ext". See \ref load_fileList
97 *
98 */
99 KLIPreduction( const std::string &dir, ///< [in] the directory to search for files.
100 const std::string &prefix, ///< [in] the prefix of the files
101 const std::string &ext ///< [in] the extension of the files
102 );
103
104 /// Construct using a file containing the target file list
105 /** Populates the \ref m_fileList vector by reading the file, which should
106 * be a single column of new-line delimited file names.
107 */
108 explicit KLIPreduction( const std::string &fileListFile /**< [in] the full path to a file
109 containing a list of files */ );
110
111 // Construct and load the target file list and the RDI file list
112 /** Populates the \ref m_fileList vector by searching on disk for files
113 * which match "dir/prefix*.ext". See \ref load_fileList
114 *
115 * Populates the \ref m_RDIfileList vector by searching on disk for files
116 * which match "RDIdir/RDIprefix*.RDIext". See \ref load_RDIfileList
117 */
118 KLIPreduction( const std::string &dir, ///< [in] the directory to search for files.
119 const std::string &prefix, ///< [in] the prefix of the files
120 const std::string &ext, ///< [in] the extension of the files, default is .fits
121 const std::string &RDIdir, ///< [in] the directory to search for the reference files.
122 const std::string &RDIprefix, /**< [in] the initial part of the file name for the
123 reference files. Can be empty "". */
124 const std::string &RDIext = "" /**< [in] [optional] the extension to append to the RDI file
125 name, must include the '.'. If empty "" then same
126 extension as target files is used. */
127 );
128
129 /// Construct using a file containing the target file list and a file
130 /// containing the RDI target file list
131 /** Populates the \ref m_fileList vector by reading the file, which should
132 * be a single column of new-line delimited file names.
133 *
134 * Populates the \ref m_RDIfileList vector by reading the file, which should
135 * be a single column of new-line delimited file names.
136 */
137 explicit KLIPreduction( const std::string &fileListFile, ///< [in] a file name path to read for
138 ///< the target file names.
139 const std::string &RDIfileListFile ///< [in] a file name path to read
140 ///< for the reference file names.
141 );
142
143 virtual ~KLIPreduction();
144
145 /// Specify how the data are centered for PCA within each search region
146 /** Can have the following values:
147 * - <b>HCI::meanSubMethod::imageMean</b> = the mean of each image (within the search
148 * region) is subtracted from itself
149 * - <b>HCI::meanSubMethod::imageMedian</b> = the median of each image (within the search
150 * region) is subtracted from itself
151 * - <b>HCI::meanSubMethod::imageMode</b> = the mode of each image (within the search
152 * region) is subtracted from itself
153 * - <b>HCI::meanSubMethod::meanImage</b> = the mean image of the data is subtracted from
154 * each image
155 * - <b>HCI::meanSubMethod::medianImage</b> = the median image of the data is subtracted
156 * from each image
157 */
159
160 /// Specify if each pixel time-series is normalized
161 /** This normalizaton is applied after centering. Can have the following values:
162 * - <b>HCI::pixelTSNormMethod::none</b>: no normalization (the default)
163 * - <b>HCI::pixelTSNormMethod::rms</b>: divide by the time-series rms
164 * - <b>HCI::pixelTSNormMethod::rmsSigmaClipped</b>: divide by the sigma-slipped time-series rms.
165 * The sigma is provided by m_pixelTSSigma.
166 */
168
169 realT m_pixelTSSigma {3}; ///< Sigma-clipping parameter for pixel time-series normalization
170
171 int m_padSize{ 4 };
172
173 /// Specifies the number of modes to include in the PSF.
174 /** The output image is a cube with a plane for each entry in m_Nmodes.
175 * Only the number of eigenvalues required for the maximum value of m_Nmodes
176 * are calculated, so this can have an important impact on speed.
177 *
178 * Can be initialized as:
179 * \code
180 * red.m_Nmodes={5,10,15,20};
181 * \endcode
182 *
183 */
184 std::vector<int> m_Nmodes;
185
186 int m_maxNmodes{ 0 };
187
188 /// Specify the minimum pixel difference at the inner edge of the search
189 /// region
190 realT m_minDPx{ 0 };
191
192 /// Specify the maximum pixel difference at the inner edge of the search
193 /// region
194 realT m_maxDPx{ 0 };
195
196 /// Controls how reference images are excluded, if at all, from the
197 /// covariance matrix for each target image based on a minimum criterion.
198 /** Can have the following values:
199 * - <b>HCI::excludeNone</b> = no exclusion, all images included [default]
200 * - <b>HCI::excludePixel</b> = exclude based on pixels of rotation at the
201 * inner edge of the region
202 * - <b>HCI::excludeAngle</b> = exclude based on degrees of rotation at the
203 * inner edge of the region
204 * - <b>HCI::excludeImno</b> = exclude based on number of images
205 */
207
208 /// Controls how reference images are excluded, if at all, from the
209 /// covariance matrix for each target image based on a maximum criterion.
210 /** Can have the following values:
211 * - <b>HCI::excludeNone</b> = no exclusion, all images included [default]
212 * - <b>HCI::excludePixel</b> = exclude based on pixels of rotation at the
213 * inner edge of the region
214 * - <b>HCI::excludeAngle</b> = exclude based on degrees of rotation at the
215 * inner edge of the region
216 * - <b>HCI::excludeImno</b> = exclude based on number of images
217 */
219
220 /// Number of reference images to include in the covariance matrix
221 /** If > 0, then at most this many images, determined by highest
222 * cross-correlation, are included. This is determined after
223 * rotational/image-number exclusion. If == 0, then all reference images are
224 * included.
225 */
227
228 /// Controls how number of included images is calculated.
229 /** The number of included images is calculated after exclusion is complete.
230 * Can have the following values:
231 * - <b>HCI::includeAll</b> = all remaining images are included [default]
232 * - <b>HCI::includeCorr</b> = the m_includeRefNum of the remaining images
233 * which are most correlated with the target are included
234 * - <b>HCI::includeTime</b> = the m_includeRefNum of the remaining images
235 * which are closest in time to the target are included
236 * - <b>HCI::includeAngle</b> = the m_includeRefNum of the remaining images
237 * which are closest in angle to the target are included
238 * - <b>HCI::includeImno</b> = the m_includeRefNum of the remaining images
239 * which are closest in image number to the target are included
240 */
242
243 eigenImage<int> m_imsIncluded;
244
245
246 bool m_rightReason = false;
247
248 realT m_rightReasonRadius = 2.5;
249
250 /// Subtract the basis mean from each of the images
251 /** The mean is subtracted according to m_meanSubMethod.
252 */
253 void meanSubtract( eigenCube<realT> &rims, ///< [in.out] The reference images. These are
254 ///< mean subtracted on output.
255 eigenCube<realT> &tims, ///< [in.out] The target images, which can be the same
256 ///< cube as rims (tested by pointer comparison), in which
257 ///< case they will be ignored. Mean subtractedon output.
258 eigenImageT &cmask, ///< [in] the cutout mask. Ignored if empty.
259 std::vector<realT> &sds ///< [out] The standard deviation of the mean
260 ///< subtracted reference images.
261 );
262
263 std::vector<_realT> m_minr;
264 std::vector<_realT> m_maxr;
265 std::vector<_realT> m_minq;
266 std::vector<_realT> m_maxq;
267
268 /// Run KLIP in a set of geometric search regions.
269 /** The arguments are 4 vectors, where each entry defines one component of
270 * the search region.
271 *
272 * \returns 0 on success
273 * \returns -1 on error
274 *
275 */
276 int regions( const std::vector<realT> & minr, ///< [in]
277 const std::vector<realT> & maxr, ///< [in]
278 const std::vector<realT> & minq, ///< [in]
279 const std::vector<realT> & maxq ///< [in]
280 );
281
282 /// Run KLIP in a geometric search region.
283 /** \overload
284 *
285 * \returns 0 on success
286 * \returns -1 on error
287 *
288 */
289 int regions( realT minr, ///< [in]
290 realT maxr, ///< [in]
291 realT minq, ///< [in]
292 realT maxq ///< [in]
293 )
294 {
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 );
299
300 return regions( vminr, vmaxr, vminq, vmaxq );
301 }
302
303 void worker( eigenCube<realT> &rims,
304 eigenCube<realT> &tims,
305 eigenImageT &cmask,
306 std::vector<size_t> &idx,
307 realT dang,
308 realT dangMax );
309
310 int finalProcess();
311
312 int processPSFSub( const std::string &dir, const std::string &prefix, const std::string &ext );
313
314 double t_worker_begin{ 0 };
315 double t_worker_end{ 0 };
316
317 double t_eigenv{ 0 };
318 double t_klim{ 0 };
319 double t_psf{ 0 };
320
321 void dump_times()
322 {
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 );
338 }
339};
340
341template <typename _realT, class _derotFunctObj, typename _evCalcT>
345
346template <typename _realT, class _derotFunctObj, typename _evCalcT>
348 const std::string &prefix,
349 const std::string &ext )
350 : ADIobservation<_realT, _derotFunctObj>( dir, prefix, ext )
351{
352}
353
354template <typename _realT, class _derotFunctObj, typename _evCalcT>
356 : ADIobservation<_realT, _derotFunctObj>( fileListFile )
357{
358}
359
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 )
368{
369}
370
371template <typename _realT, class _derotFunctObj, typename _evCalcT>
373 const std::string &RDIfileListFile )
374 : ADIobservation<_realT, _derotFunctObj>( fileListFile, RDIfileListFile )
375{
376}
377
378template <typename _realT, class _derotFunctObj, typename _evCalcT>
380{
381}
382
383template <typename realT, class derotFunctObj, typename evCalcT>
385 eigenCube<realT> &tims,
386 eigenImageT &cmask,
387 std::vector<realT> &norms )
388{
389
390 norms.resize( rims.planes() );
391
392 bool haveMask = false;
393 if( cmask.rows() > 0 && cmask.cols() > 0 )
394 {
395 haveMask = true;
396 }
397
398 if( m_meanSubMethod == HCI::meanSubMethod::meanImage || m_meanSubMethod == HCI::meanSubMethod::medianImage )
399 {
400 eigenImageT mean;
401
402 if( m_meanSubMethod == HCI::meanSubMethod::meanImage )
403 {
404 rims.mean( mean );
405 }
406 else if( m_meanSubMethod == HCI::meanSubMethod::medianImage )
407 {
408 rims.median( mean );
409 }
410
411 for( int n = 0; n < rims.planes(); ++n )
412 {
413 rims.image( n ) -= mean;
414
415 if( haveMask )
416 {
417 rims.image( n ) *= cmask;
418 }
419
420 realT immean = rims.image( n ).mean();
421 norms[n] = ( rims.image( n ) - immean ).matrix().norm();
422 }
423
424 if( &tims != &rims )
425 {
426 if( m_meanSubMethod == HCI::meanSubMethod::meanImage )
427 {
428 tims.mean( mean );
429 }
430 else if( m_meanSubMethod == HCI::meanSubMethod::medianImage )
431 {
432 tims.median( mean );
433 }
434
435 for( int n = 0; n < tims.planes(); ++n )
436 {
437 tims.image( n ) -= mean;
438
439 if( haveMask )
440 {
441 tims.image( n ) *= cmask;
442 }
443 }
444 }
445 }
446 else
447 {
448 realT mean;
449 std::vector<realT> work; // Working memmory for median calc
450
451 for( int n = 0; n < rims.planes(); ++n )
452 {
453 if( m_meanSubMethod == HCI::meanSubMethod::imageMean )
454 {
455 mean = rims.image( n ).mean();
456 }
457 else if( m_meanSubMethod == HCI::meanSubMethod::imageMedian )
458 {
459 mean = imageMedian( rims.image( n ), &work );
460 }
461
462 rims.image( n ) -= mean;
463
464 if( haveMask )
465 {
466 rims.image( n ) *= cmask;
467 }
468 // Because we might not have used the mean, we need to re-mean to
469 // make this the standard deviation
470 realT immean = rims.image( n ).mean();
471 norms[n] = ( rims.image( n ) - immean ).matrix().norm();
472 }
473
474 if( &tims != &rims )
475 {
476 for( int n = 0; n < tims.planes(); ++n )
477 {
478 if( m_meanSubMethod == HCI::meanSubMethod::imageMean )
479 {
480 mean = tims.image( n ).mean();
481 }
482 else if( m_meanSubMethod == HCI::meanSubMethod::imageMedian )
483 {
484 mean = imageMedian( tims.image( n ), &work );
485 }
486
487 tims.image( n ) -= mean;
488 if( haveMask )
489 {
490 tims.image( n ) *= cmask;
491 }
492 }
493 }
494 }
495
496 if(m_pixelTSNormMethod == HCI::pixelTSNormMethod::unknown)
497 {
498 mxThrowException(err::invalidconfig, "KlipReduction::meanSubtract", "pixelTSNormMethod is unknown");
499 }
500 if(m_pixelTSNormMethod == HCI::pixelTSNormMethod::rmsSigmaClipped)
501 {
502 mxThrowException(err::invalidconfig, "KlipReduction::meanSubtract", "pixelTSNormMethod is rmsSigmaClipped, which is not implemented");
503 }
504
505 if( m_pixelTSNormMethod != HCI::pixelTSNormMethod::none )
506 {
507 std::cerr << "normalizing pixels\n";
508 std::vector<realT> pixs(rims.planes());
509
510 for(int cc = 0; cc < rims.cols(); ++cc)
511 {
512 for(int rr = 0; rr < rims.rows(); ++rr)
513 {
514 if(haveMask)
515 {
516 if(cmask(rr,cc) == 0)
517 {
518 continue;
519 }
520 }
521
522 //We bother to load a vector in prep to add sigma clipping later.
523 for(int pp = 0; pp < rims.planes(); ++pp)
524 {
525 pixs[pp] = rims.image(pp)(rr,cc);
526 }
527
528 realT sd = sqrt(math::vectorVariance(pixs));
529
530 for(int pp = 0; pp < rims.planes(); ++pp)
531 {
532 rims.image(pp)(rr,cc) /= sd;
533 }
534 }
535 }
536 }
537}
538
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 )
544{
545 this->t_begin = sys::get_curr_time();
546
547 m_minr = minr;
548 m_maxr = maxr;
549 m_minq = minq;
550 m_maxq = maxq;
551
552 m_maxNmodes = m_Nmodes[0];
553 for( size_t i = 1; i < m_Nmodes.size(); ++i )
554 {
555 if( m_Nmodes[i] > m_maxNmodes )
556 m_maxNmodes = m_Nmodes[i];
557 }
558
559 std::cerr << "Beginning\n";
560
561 if( this->m_imSize == 0 )
562 {
563 this->m_imSize = 2 * ( *std::max_element( m_maxr.begin(), m_maxr.end() ) + m_padSize );
564
565 std::cerr << "set image size based on regions to " << this->m_imSize << "\n";
566 }
567
568 if( !this->m_filesRead )
569 {
570 if( this->readFiles() < 0 )
571 return -1;
572 }
573
574 // CHECK IF RDI HERE
575 if( !this->m_RDIfilesRead && this->m_RDIfileList.size() != 0 )
576 {
577 if( this->readRDIFiles() < 0 )
578 {
579 return -1;
580 }
581 }
582
583 if( this->m_preProcess_only && !this->m_skipPreProcess )
584 {
585 std::cerr << "Pre-processing complete, stopping.\n";
586 return 0;
587 }
588
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 )
592 {
593 this->m_psfsub[n].resize( this->m_Nrows, this->m_Ncols, this->m_Nims );
594 this->m_psfsub[n].cube().setZero();
595 }
596
597 // Make radius and angle images
598 eigenImageT rIm( this->m_Nrows, this->m_Ncols );
599 eigenImageT qIm( this->m_Nrows, this->m_Ncols );
600
601 radAngImage<math::degreesT<realT>>( rIm, qIm, .5 * ( this->m_Nrows - 1 ), .5 * ( this->m_Ncols - 1 ) );
602
603 m_imsIncluded.resize( this->m_Nims, this->m_Nims );
604 m_imsIncluded.setConstant( 1 );
605
606 if( this->m_refIms.planes() > 0 ) // RDI
607 {
608 std::cerr << "******* RDI MODE **********\n";
609 }
610 else // ADI
611 {
612 std::cerr << "******* ADI MODE **********\n";
613 }
614
615 std::cerr << "processing " << minr.size() << " regions\n";
616
617 //******** For each region do this:
618 for( size_t regno = 0; regno < minr.size(); ++regno )
619 {
620 std::cerr << " region " << regno+1 << ": " << m_minr[regno] << "-" << m_maxr[regno] << " pixels, ";
621 std::cerr << m_minq[regno] << "-" << m_maxq[regno] << " degrees. \n";
622
623 eigenImageT *maskPtr = nullptr;
624
625 if( this->m_mask.rows() == this->m_Nrows && this->m_mask.cols() == this->m_Ncols )
626 {
627 maskPtr = &this->m_mask;
628 }
629
630 std::vector<size_t> idx = annulusIndices<math::degreesT<realT>>( rIm,
631 qIm,
632 .5 * ( this->m_Nrows - 1 ),
633 .5 * ( this->m_Ncols - 1 ),
634 m_minr[regno],
635 m_maxr[regno],
636 m_minq[regno],
637 m_maxq[regno],
638 maskPtr );
639
640 // Create storage for the R-ims and psf-subbed Ims
641 eigenCube<realT> tims( idx.size(), 1, this->m_Nims );
642
643 //------If doing RDI, create bims
644 eigenCube<realT> rims;
645
646 //------Get the mask cutout too
647 eigenImageT cmask;
648
649 if( this->m_refIms.planes() > 0 )
650 {
651 rims.resize( idx.size(), 1, this->m_Nims );
652 }
653
654 for( int i = 0; i < this->m_Nims; ++i )
655 {
656 auto tim = tims.image( i );
657 cutImageRegion( tim, this->m_tgtIms.image( i ), idx, false );
658 }
659
660 for( int p = 0; p < this->m_refIms.planes(); ++p )
661 {
662 auto rim = rims.image( p );
663 cutImageRegion( rim, this->m_refIms.image( p ), idx, false );
664 }
665
666 if( this->m_maskFile != "" )
667 {
668 cutImageRegion( cmask, this->m_mask, idx, true );
669 }
670
671 #if 0
672
673 std::vector<realT> sds;
674
675 //*** First mean subtract ***//
676 meanSubtract( tims, tims, cmask, sds );
677
678
680 ffff.write("tims.fits", tims);
681
682 std::ofstream fout("idx.dat");
683 std::ofstream aout("derot.dat");
684 for(auto i : idx)
685 {
686 fout << i << "\n";
687 }
688 fout.close();
689
690 for(int pp =0; pp < tims.planes(); ++pp)
691 {
692 aout << this->m_derotF.derotAngle( pp ) << "\n";
693 }
694 aout.close();
695
696 fout.open("dims.dat");
697 fout << this->m_Nrows << "\n";
698 fout << this->m_Ncols << "\n";
699 fout.close();
700
701
702 exit(0);
703 #endif
704
705 realT dang = 0;
706 realT dangMax = 0;
707
708 if( m_minDPx < 0 )
709 {
710 m_excludeMethod = HCI::excludeNone;
711 }
712 if( m_maxDPx < 0 )
713 {
714 m_excludeMethodMax = HCI::excludeNone;
715 }
716
717 //------- If doing RDI, excludeMethod and excludeMethodMax must be none!
718 if( this->m_refIms.planes() > 0 )
719 {
720 m_excludeMethod = HCI::excludeNone;
721 m_excludeMethodMax = HCI::excludeNone;
722 }
723
724 if( m_excludeMethod == HCI::excludePixel )
725 {
726 dang = fabs( atan( m_minDPx / minr[regno] ) );
727 }
728 else if( m_excludeMethod == HCI::excludeAngle )
729 {
730 dang = math::dtor( m_minDPx );
731 }
732 else if( m_excludeMethod == HCI::excludeImno )
733 {
734 dang = m_minDPx;
735 }
736
737 if( m_excludeMethodMax == HCI::excludePixel )
738 {
739 dangMax = fabs( atan( m_maxDPx / minr[regno] ) );
740 }
741 else if( m_excludeMethodMax == HCI::excludeAngle )
742 {
743 dangMax = math::dtor( m_maxDPx );
744 }
745 else if( m_excludeMethodMax == HCI::excludeImno )
746 {
747 dangMax = m_maxDPx;
748 }
749
750 //------- If doing RDI, call this with rims, bims
751 //*** Dispatch the work
752 if( this->m_refIms.planes() > 0 ) // RDI
753 {
754 worker( rims, tims, cmask, idx, dang, dangMax );
755 }
756 else // ADI
757 {
758 worker( tims, tims, cmask, idx, dang, dangMax );
759 }
760 }
761
763 ffii.write( "imsIncluded.fits", m_imsIncluded );
764
765 if( finalProcess() < 0 )
766 {
767 std::cerr << "Error in final processing\n";
768 }
769
770 this->t_end = sys::get_curr_time();
771
772 dump_times();
773
774 return 0;
775}
776
777struct cvEntry
778{
779 int index;
780 double cvVal;
781 double angle;
782 bool included{ true };
783};
784
785template <typename eigenT, typename eigenTin>
786void extractRowsAndCols( eigenT &out, const eigenTin &in, const std::vector<size_t> &idx )
787{
788 out.resize( idx.size(), idx.size() );
789
790 for( size_t i = 0; i < idx.size(); ++i )
791 {
792 for( size_t j = 0; j < idx.size(); ++j )
793 {
794 out( i, j ) = in( idx[i], idx[j] );
795 }
796 }
797}
798
799template <typename eigenT, typename eigenTin>
800void extractCols( eigenT &out, const eigenTin &in, const std::vector<size_t> &idx )
801{
802 out.resize( in.rows(), idx.size() );
803
804 for( size_t i = 0; i < idx.size(); ++i )
805 {
806 out.col( i ) = in.col( idx[i] ); // it1->index);
807 }
808}
809
810template <typename realT, typename eigenT, typename eigenTv, class derotFunctObj>
811void collapseCovar( eigenT &cutCV,
812 const eigenT &CV,
813 const std::vector<realT> &sds,
814 eigenT &rimsCut,
815 const eigenTv &rims,
816 int imno,
817 double dang,
818 double dangMax,
819 int Nims,
820 int excludeMethod,
821 int excludeMethodMax,
822 int includeRefNum,
823 const derotFunctObj &derotF,
824 eigenImage<int> &imsIncluded )
825{
826 std::vector<cvEntry> allidx( Nims );
827
828 // Initialize the vector of cvEntries
829 for( int i = 0; i < Nims; ++i )
830 {
831 allidx[i].index = i;
832 allidx[i].angle = derotF.derotAngle( i );
833
834 // CV is lower-triangular
835 if( i <= imno )
836 {
837 allidx[i].cvVal = CV( imno, i ) / ( sds[i] * sds[imno] );
838 }
839 else
840 {
841 allidx[i].cvVal = CV( i, imno ) / ( sds[i] * sds[imno] );
842 }
843 }
844
845 if( excludeMethod == HCI::excludePixel || excludeMethod == HCI::excludeAngle )
846 {
847 for( size_t j = 0; j < Nims; ++j )
848 {
849 if( fabs( math::angleDiff<math::radiansT<realT>>( derotF.derotAngle( j ), derotF.derotAngle( imno ) ) ) <=
850 dang )
851 {
852 allidx[j].included = false;
853 }
854 }
855 }
856 else if( excludeMethod == HCI::excludeImno )
857 {
858 for( size_t j = 0; j < Nims; ++j )
859 {
860 if( fabs( (long)j - imno ) <= dang )
861 {
862 allidx[j].included = false;
863 }
864 }
865 }
866
867 if( excludeMethodMax == HCI::excludePixel || excludeMethodMax == HCI::excludeAngle )
868 {
869 for( size_t j = 0; j < Nims; ++j )
870 {
871 if( fabs( math::angleDiff<math::radiansT<realT>>( derotF.derotAngle( j ), derotF.derotAngle( imno ) ) ) >
872 dangMax )
873 allidx[j].included = false;
874 }
875 }
876 else if( excludeMethodMax == HCI::excludeImno )
877 {
878 for( size_t j = 0; j < Nims; ++j )
879 {
880 if( fabs( (long)j - imno ) > dangMax )
881 allidx[j].included = false;
882 }
883 }
884
885 if( includeRefNum > 0 && (size_t)includeRefNum < allidx.size() )
886 {
887 long kept = 0;
888 for( size_t j = 0; j < Nims; ++j )
889 {
890 if( allidx[j].included == true )
891 {
892 ++kept;
893 }
894 }
895
896 // Get a vector for sorting
897 std::vector<realT> cvVal;
898 cvVal.resize( kept );
899 size_t k = 0;
900 for( size_t j = 0; j < Nims; ++j )
901 {
902 if( allidx[j].included == true )
903 {
904 cvVal[k] = allidx[j].cvVal;
905 ++k;
906 }
907 }
908
909 // Partially sort the correlation values
910 std::nth_element( cvVal.begin(), cvVal.begin() + ( kept - includeRefNum ), cvVal.end() );
911
912 realT mincorr = cvVal[kept - includeRefNum];
913 //std::cerr << " Minimum correlation: " << mincorr << "\n";
914
915 for( size_t j = 0; j < Nims; ++j )
916 {
917 if( allidx[j].cvVal < mincorr )
918 {
919 allidx[j].included = false;
920 }
921 }
922 }
923
924 std::vector<size_t> keepidx;
925 for( size_t j = 0; j < Nims; ++j )
926 {
927 imsIncluded( imno, j ) = allidx[j].included;
928
929 if( allidx[j].included )
930 keepidx.push_back( j );
931 }
932
933 // std::cerr << " Keeping " << keepidx.size() << " reference images out of
934 // "
935 // << Nims << " (" << Nims-keepidx.size() << " rejected)\n";
936
937 if( keepidx.size() == 0 )
938 {
939 std::cerr << "\n\n" << imno << "\n\n";
940 }
941
942 extractRowsAndCols( cutCV, CV, keepidx );
943 extractCols( rimsCut, rims, keepidx );
944}
945
946template <typename _realT, class _derotFunctObj, typename _evCalcT>
948 eigenCube<_realT> &tims,
949 eigenImageT &cmask,
950 std::vector<size_t> &idx,
951 realT dang,
952 realT dangMax )
953{
954
955 t_worker_begin = sys::get_curr_time();
956
957 std::vector<realT> sds;
958
959 eigenImageT meanim;
960
961 //*** First mean subtract ***//
962 meanSubtract( rims, tims, cmask, sds );
963
964 //*** Form lower-triangle covariance matrix
965 eigenImageT cv;
966
967 math::eigenSYRK( cv, rims.cube() );
968
970 ff.write( "cv.fits", cv );
971 ipc::ompLoopWatcher<> status( this->m_Nims, std::cerr );
972
973 // Pre-calculate KL images once if we are exclude none OR IF RDI
974 eigenImageT master_klims;
975 eigenImageT master_projMat;
976
977 eigenImage<realT> rrMask;
978 if( m_rightReason )
979 {
980 rrMask.resize( tims.rows(), tims.rows() );
981 rrMask.setConstant( 1 );
982
983 eigenImageT rrmim;
984 rrmim.resize( this->m_Nrows, this->m_Ncols );
985 for( int rr = 0; rr < rrMask.rows(); ++rr )
986 {
987 rrmim.setConstant( 1 );
988
989 //Calculate the 2D coords of this pixel
990 int jj = idx[rr] / this->m_Ncols;
991 int ii = idx[rr] - this->m_Ncols * jj;
992
993 maskCircle( rrmim, ii, jj, m_rightReasonRadius, 0, 0 );
994
995 //Extract the 2D r.r. mask into the row-vector image.
996 for( int cc = 0; cc < rrMask.cols(); ++cc )
997 {
998 rrMask( rr, cc ) = rrmim.data()[idx[cc]];
999 }
1000 }
1001
1002 ff.write( "rrMask.fits", rrMask );
1003 }
1004
1005 if( m_excludeMethod == HCI::excludeNone && m_excludeMethodMax == HCI::excludeNone && m_includeRefNum == 0 )
1006 {
1007 double teigenv;
1008 double tklim;
1009
1010 math::calcKLModes<double>( master_klims, cv, rims.cube(), m_maxNmodes, nullptr, &teigenv, &tklim );
1011
1012 if( m_rightReason )
1013 {
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 );
1018 }
1019
1020 t_eigenv += teigenv;
1021 t_klim += tklim;
1022 }
1023
1024 // clang-format off
1025 #pragma omp parallel // num_threads(20)
1026 //clang-format on
1027 {
1028 // We need local copies for each thread. Only way this works, for
1029 // whatever reason.
1030 eigenImageT cfs; // The coefficients
1031 eigenImageT psf;
1032 eigenImageT rims_cut;
1033 eigenImageT cv_cut;
1034 eigenImageT klims;
1035 eigenImageT projMat;
1036
1038
1039 if( m_excludeMethod == HCI::excludeNone && m_excludeMethodMax == HCI::excludeNone &&
1040 m_includeRefNum == 0 ) // OR RDI
1041 {
1042 klims = master_klims;
1043 if( m_rightReason )
1044 {
1045 projMat = master_projMat;
1046 }
1047 }
1048
1049 // clang-format off
1050 #pragma omp for
1051 // clang-format on
1052 for( int imno = 0; imno < this->m_Nims; ++imno )
1053 {
1054 if( m_excludeMethod != HCI::excludeNone || m_excludeMethodMax != HCI::excludeNone || m_includeRefNum != 0 )
1055 {
1056 collapseCovar<realT>( cv_cut,
1057 cv,
1058 sds,
1059 rims_cut,
1060 rims.asVectors(),
1061 imno,
1062 dang,
1063 dangMax,
1064 this->m_Nims,
1065 this->m_excludeMethod,
1066 this->m_excludeMethodMax,
1067 this->m_includeRefNum,
1068 this->m_derotF,
1069 m_imsIncluded );
1070
1071 /**** Now calculate the K-L Images ****/
1072 double teigenv, tklim;
1073 math::calcKLModes( klims, cv_cut, rims_cut, m_maxNmodes, &mem, &teigenv, &tklim );
1074
1075 if( m_rightReason )
1076 {
1077 projMat = ( klims.matrix().transpose() * klims.matrix() ).array();
1078 projMat *= rrMask;
1079 }
1080
1081 t_eigenv += teigenv;
1082 t_klim += tklim;
1083 }
1084 cfs.resize( 1, klims.rows() );
1085
1086 double t0 = sys::get_curr_time();
1087
1088 if( !m_rightReason )
1089 {
1090 for( int j = 0; j < cfs.size(); ++j )
1091 {
1092 cfs( j ) = klims.row( j ).matrix().dot( tims.cube().col( imno ).matrix() );
1093 }
1094
1095 for( size_t mode_i = 0; mode_i < m_Nmodes.size(); ++mode_i )
1096 {
1097 psf = cfs( cfs.size() - 1 ) * klims.row( cfs.size() - 1 );
1098
1099 // Count down, since eigenvalues are returned in increasing
1100 // order
1101 // handle case where cfs.size() < m_Nmodes[mode_i], i.e.
1102 // when more modes than images.
1103 for( int j = cfs.size() - 2; j >= cfs.size() - m_Nmodes[mode_i] && j >= 0; --j )
1104 {
1105 psf += cfs( j ) * klims.row( j );
1106 }
1107
1108 // #pragma omp critical
1110 this->m_psfsub[mode_i].cube().col( imno ), tims.cube().col( imno ) - psf.transpose(), idx );
1111 }
1112 }
1113 else
1114 {
1115 psf = projMat.matrix() * tims.cube().col( imno ).matrix();
1116 insertImageRegion( this->m_psfsub[0].cube().col( imno ), tims.cube().col( imno ) - psf, idx );
1117 }
1118
1119 t_psf += ( sys::get_curr_time() - t0 ); /// omp_get_num_threads();
1120
1121 status.incrementAndOutputStatus();
1122
1123 } // for imno
1124 } // openmp parrallel
1125
1126 t_worker_end = sys::get_curr_time();
1127}
1128
1129template <typename _realT, class _derotFunctObj, typename _evCalcT>
1131{
1132 if( this->m_postMedSub )
1133 {
1134 std::cerr << "Subtracting medians in post\n";
1135
1136 for( size_t n = 0; n < this->m_psfsub.size(); ++n )
1137 {
1138 // clang-format off
1139 #pragma omp parallel
1140 // clang-format on
1141 {
1142 eigenImage<realT> medim;
1143
1144 this->m_psfsub[n].median( medim );
1145
1146 // clang-format off
1147 #pragma omp for
1148 // clang-format on
1149 for( int i = 0; i < this->m_psfsub[n].planes(); ++i )
1150 {
1151 this->m_psfsub[n].image( i ) -= medim;
1152 }
1153 }
1154 }
1155 }
1156
1157 if( this->m_doDerotate )
1158 {
1159 std::cerr << "derotating\n";
1160 this->derotate();
1161 }
1162
1163 if( this->m_combineMethod > 0 )
1164 {
1165 std::cerr << "combining\n";
1166 this->combineFinim();
1167 }
1168
1169 if( this->m_doWriteFinim == true || this->m_doOutputPSFSub == true )
1170 {
1171 std::cerr << "writing\n";
1172
1173 fits::fitsHeader head;
1174
1175 this->ADIobservation<_realT, _derotFunctObj>::stdFitsHeader( &head );
1176
1177 head.append( "", fits::fitsCommentType(), "----------------------------------------" );
1178 head.append( "", fits::fitsCommentType(), "mx::KLIPreduction parameters:" );
1179 head.append( "", fits::fitsCommentType(), "----------------------------------------" );
1180
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" );
1183
1184 std::stringstream str;
1185
1186 if( m_Nmodes.size() > 0 )
1187 {
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" );
1192 }
1193
1194 head.append<bool>( "RIGHT REASON", m_rightReason, "whether or not the right reason mask is used");
1195 if( m_rightReason)
1196 {
1197 head.append<realT>( "RIGHT REASON RADIUS", m_rightReasonRadius, "radius of the right reason mask");
1198 }
1199
1200 if( m_minr.size() > 0 )
1201 {
1202 str.str( "" );
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)" );
1207 }
1208
1209 if( m_maxr.size() > 0 )
1210 {
1211 str.str( "" );
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)" );
1216 }
1217
1218 if( m_minq.size() > 0 )
1219 {
1220 str.str( "" );
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)" );
1225 }
1226
1227 if( m_maxq.size() > 0 )
1228 {
1229 str.str( "" );
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)" );
1234 }
1235
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)" );
1238
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)" );
1241
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" );
1244
1245 if( this->m_doWriteFinim == true && this->m_combineMethod > 0 )
1246 {
1247 this->writeFinim( &head );
1248 }
1249
1250 if( this->m_doOutputPSFSub )
1251 {
1252 this->outputPSFSub( &head );
1253 }
1254 }
1255
1256 return 0;
1257}
1258
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 )
1263
1264{
1265 std::cerr << "Beginning PSF Subtracted Image Processing\n";
1266
1267 // Load first file to condigure based on its header.
1268 std::vector<std::string> flist = ioutils::getFileNames( dir, prefix, "000", ext );
1269
1271 eigenImage<realT> im;
1273
1274 ff.read( im, fh, flist[0] );
1275
1276 if( fh.count( "MEANSUBM" ) == 0 )
1277 {
1278 mxError( "KLIPReduction", MXE_PARAMNOTSET, "MEANSUBM not found in FITS header." );
1279 return -1;
1280 }
1281 m_meanSubMethod = HCI::meanSubMethodFmStr( fh["MEANSUBM"].String() );
1282 std::cerr << "meanSubMethod: " << HCI::meanSubMethodStr( m_meanSubMethod ) << "\n";
1283
1284 if( fh.count( "NMODES" ) == 0 )
1285 {
1286 mxError( "KLIPReduction", MXE_PARAMNOTSET, "NMODES not found in FITS header." );
1287 return -1;
1288 }
1289 ioutils::parseStringVector( m_Nmodes, fh["NMODES"].String(), "," );
1290 if( m_Nmodes.size() == 0 )
1291 {
1292 mxError( "KLIPReduction", MXE_PARSEERR, "NMODES vector did not parse correctly." );
1293 return -1;
1294 }
1295 std::cerr << "nModes: " << fh["NMODES"].String() << "\n";
1296
1297 /* -- REGMINR -- */
1298 if( fh.count( "REGMINR" ) == 0 )
1299 {
1300 mxError( "KLIPReduction", MXE_PARAMNOTSET, "REGMINR not found in FITS header." );
1301 return -1;
1302 }
1303 ioutils::parseStringVector( m_minr, fh["REGMINR"].String(), "," );
1304 if( m_minr.size() == 0 )
1305 {
1306 mxError( "KLIPReduction", MXE_PARSEERR, "REGMINR vector did not parse correctly." );
1307 return -1;
1308 }
1309 std::cerr << "minr: " << fh["REGMINR"].String() << "\n";
1310
1311 /* -- REGMAXR -- */
1312 if( fh.count( "REGMAXR" ) == 0 )
1313 {
1314 mxError( "KLIPReduction", MXE_PARAMNOTSET, "REGMAXR not found in FITS header." );
1315 return -1;
1316 }
1317 ioutils::parseStringVector( m_maxr, fh["REGMAXR"].String(), "," );
1318 if( m_maxr.size() == 0 )
1319 {
1320 mxError( "KLIPReduction", MXE_PARSEERR, "REGMAXR vector did not parse correctly." );
1321 return -1;
1322 }
1323 std::cerr << "minr: " << fh["REGMAXR"].String() << "\n";
1324
1325 /* -- REGMINQ -- */
1326 if( fh.count( "REGMINQ" ) == 0 )
1327 {
1328 mxError( "KLIPReduction", MXE_PARAMNOTSET, "REGMINQ not found in FITS header." );
1329 return -1;
1330 }
1331 ioutils::parseStringVector( m_minq, fh["REGMINQ"].String(), "," );
1332 if( m_minq.size() == 0 )
1333 {
1334 mxError( "KLIPReduction", MXE_PARSEERR, "REGMINQ vector did not parse correctly." );
1335 return -1;
1336 }
1337 std::cerr << "minq: " << fh["REGMINQ"].String() << "\n";
1338
1339 /* -- REGMAXR -- */
1340 if( fh.count( "REGMAXR" ) == 0 )
1341 {
1342 mxError( "KLIPReduction", MXE_PARAMNOTSET, "REGMAXR not found in FITS header." );
1343 return -1;
1344 }
1345 ioutils::parseStringVector( m_maxq, fh["REGMAXR"].String(), "," );
1346 if( m_maxq.size() == 0 )
1347 {
1348 mxError( "KLIPReduction", MXE_PARSEERR, "REGMAXR vector did not parse correctly." );
1349 return -1;
1350 }
1351 std::cerr << "minr: " << fh["REGMAXR"].String() << "\n";
1352
1353 if( fh.count( "EXMTHDMN" ) == 0 )
1354 {
1355 mxError( "KLIPReduction", MXE_PARAMNOTSET, "EXMTHDMN not found in FITS header." );
1356 return -1;
1357 }
1358 m_excludeMethod = HCI::excludeMethodFmStr( fh["EXMTHDMN"].String() );
1359 std::cerr << "excludeMethod: " << HCI::excludeMethodStr( m_excludeMethod ) << "\n";
1360
1361 if( fh.count( "MINDPX" ) == 0 )
1362 {
1363 mxError( "KLIPReduction", MXE_PARAMNOTSET, "MINDPX not found in FITS header." );
1364 return -1;
1365 }
1366 m_minDPx = fh["MINDPX"].value<realT>();
1367 std::cerr << "minDPx: " << m_minDPx << "\n";
1368
1369 if( fh.count( "EXMTHDMX" ) == 0 )
1370 {
1371 mxError( "KLIPReduction", MXE_PARAMNOTSET, "EXMTHDMX not found in FITS header." );
1372 return -1;
1373 }
1374 m_excludeMethodMax = HCI::excludeMethodFmStr( fh["EXMTHDMX"].String() );
1375 std::cerr << "excludeMethodMax: " << HCI::excludeMethodStr( m_excludeMethodMax ) << "\n";
1376
1377 if( fh.count( "MAXDPX" ) == 0 )
1378 {
1379 mxError( "KLIPReduction", MXE_PARAMNOTSET, "MAXDPX not found in FITS header." );
1380 return -1;
1381 }
1382 m_maxDPx = fh["MAXDPX"].value<realT>();
1383 std::cerr << "maxDPx: " << m_maxDPx << "\n";
1384
1385 if( fh.count( "INMTHDMX" ) == 0 )
1386 {
1387 mxError( "KLIPReduction", MXE_PARAMNOTSET, "INMTHDMX not found in FITS header." );
1388 return -1;
1389 }
1390 m_includeMethod = HCI::includeMethodFmStr( fh["INMTHDMX"].String() );
1391 std::cerr << "includeMethod: " << HCI::includeMethodStr( m_includeMethod ) << "\n";
1392
1393 if( fh.count( "INCLREFN" ) == 0 )
1394 {
1395 mxError( "KLIPReduction", MXE_PARAMNOTSET, "INCLREFN not found in FITS header." );
1396 return -1;
1397 }
1398 m_includeRefNum = fh["INCLREFN"].value<int>();
1399 std::cerr << "includedRefNum: " << m_includeRefNum << "\n";
1400
1401 this->m_skipPreProcess = true;
1402
1403 this->m_keywords.clear();
1404
1405 this->readPSFSub( dir, prefix, ext, m_Nmodes.size() );
1406
1407 finalProcess();
1408
1409 return 0;
1410}
1411
1412///@}
1413
1414template <typename realT>
1415class ADIDerotator;
1416
1417extern template struct KLIPreduction<float, ADIDerotator<float>, float>;
1418extern template struct KLIPreduction<float, ADIDerotator<float>, double>;
1419extern template struct KLIPreduction<double, ADIDerotator<double>, double>;
1420
1421} // namespace improc
1422} // namespace mx
1423
1424#endif //__KLIPreduction_hpp__
Defines the ADI high contrast imaging data type.
mxException for invalid config settings
Class to manage interactions with a FITS file.
Definition fitsFile.hpp:52
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.
Definition fitsFile.hpp:829
Class to manage a complete fits header.
size_t count(const std::string &keyword)
Get number of cards with a given keyword.
void clear()
Clear all cards from the header.
void append(const fitsHeaderCard &card)
Append a fitsHeaderCard to the end of the header.
An image cube with an Eigen-like API.
Definition eigenCube.hpp:30
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.
Definition constants.hpp:69
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.
Definition geo.hpp:197
realT dtor(realT q)
Convert from degrees to radians.
Definition geo.hpp:138
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.
Definition timeUtils.hpp:90
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
The mxlib c++ namespace.
Definition mxError.hpp:106
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.
int m_includeMethod
Controls how number of included images is calculated.
Type holding constants related to angle calculations in degrees.
Definition geo.hpp:87
A struct to hold the working memory for eigenSYEVR and maintain it between calls if desired.