mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
Loading...
Searching...
No Matches
HCIobservation.hpp
Go to the documentation of this file.
1/** \file HCIobservation.hpp
2 * \author Jared R. Males
3 * \brief Defines the basic high contrast imaging data type.
4 * \ingroup hc_imaging_files
5 * \ingroup image_processing_files
6 *
7 */
8
9#ifndef __HCIobservation_hpp__
10#define __HCIobservation_hpp__
11
12#include <vector>
13#include <map>
14#include <string>
15#include <fstream>
16
17#include <sys/stat.h>
18
19#include "../mxlib.hpp"
20
21#include "../mxException.hpp"
22
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"
29
30#include "eigenImage.hpp"
31#include "eigenCube.hpp"
32#include "imageFilters.hpp"
33#include "imageMasks.hpp"
34#include "imageTransforms.hpp"
35#include "imageUtils.hpp"
36
37namespace mx
38{
39
40namespace improc
41{
42
43/// Namespace for high contrast imaging enums.
44/** \ingroup hc_imaging_enums
45 */
46namespace HCI
47{
48
49/// Mean subtraction methods
50/** These control how the data in each search region is centered to meet the PCA
51 * requirement. \ingroup hc_imaging_enums
52 */
53enum class meanSubMethod
54{
55 none, ///< No mean subtraction
56 meanImage, ///< The mean image of the data is subtracted from each image
57 medianImage, ///< The median image of the data is subtracted from each image
58 imageMean, ///< The mean of each image (within the search region) is
59 ///< subtracted from itself
60 imageMedian, ///< The median of each image (within the search region) is
61 ///< subtracted from itself
62 imageMode, ///< The mode of each image (within the search region) is
63 ///< subtracted from itself
64 unknown = -1 ///< unknown value, an error
65};
66
67std::string meanSubMethodStr( meanSubMethod method );
68
69meanSubMethod meanSubMethodFmStr( const std::string &method );
70
72{
73 none, ///< no pixel time series norm
74 rms, ///< the rms of the pixel time series
75 rmsSigmaClipped, ///< the sigma clipped rms of the pixel time series
76 unknown = -1 ///< unknown value, an error
77};
78
79std::string pixelTSNormMethodStr( pixelTSNormMethod method );
80
81pixelTSNormMethod pixelTSNormMethodFmStr( const std::string &method );
82
83/// Possible combination methods
84/** \ingroup hc_imaging_enums
85 */
87{
88 noCombine, ///< Do not combine the images.
89 medianCombine, ///< Combine with the median.
90 meanCombine, ///< Combine with the mean.
91 sigmaMeanCombine, ///< Combine with the sigma clipped mean.
92 debug
93};
94
95/// Get the string name of the combineMethod integer
96/**
97 * \returns a string with the name of the combineMethod
98 */
99std::string combineMethodStr( int method /**< [in] one of the combineMethods enum members */ );
100
101/// Get the combineMethod from the corresponding string name
102/**
103 * \returns the combineMethods enum member corresponding to the string name.
104 */
105int combineMethodFmStr( const std::string &method /**< [in] the string name of the combineMethod */ );
106
107} // namespace HCI
108
109/// The basic high contrast imaging data type
110/** This class manages file reading, resizing, co-adding, pre-processing (masking and filtering),
111 * and final image combination.
112 *
113 * \tparam _realT is the floating point type in which to do all arithmetic.
114 *
115 * \ingroup hc_imaging
116 */
117template <typename _realT>
119{
120
121 /// The arithmetic type used for calculations. Does not have to match the type in images on disk.
122 typedef _realT realT;
123
124 /// The Eigen image array type basted on realT
125 typedef Eigen::Array<realT, Eigen::Dynamic, Eigen::Dynamic> eigenImageT;
126
127 ///\name Construction and Initialization
128 /** @{
129 */
130 /// Default c'tor
132
133 /// Construct and load the target file list.
134 /** Populates the \ref m_fileList vector by searching on disk for files which match
135 * "dir/prefix*.ext". See \ref load_fileList
136 *
137 */
138 HCIobservation( const std::string &dir, ///< [in] the directory to search.
139 const std::string &prefix, ///< [in] the initial part of the file name. Can be empty "".
140 const std::string &ext ///< [in] the extension to append to the file name, must include the '.'.
141 );
142
143 /// Construct using a file containing the target file list
144 /** Populates the \ref m_fileList vector by reading the file, which should be a single
145 * column of new-line delimited file names.
146 */
147 explicit HCIobservation( const std::string &fileListFile /**< [in] a file name path to read.*/ );
148
149 /// Construct and load the target file list and the RDI file list
150 /** Populates the \ref m_fileList vector by searching on disk for files which match
151 * "dir/prefix*.ext". See \ref load_fileList
152 *
153 * Populates the \ref m_RDIfileList vector by searching on disk for files which match
154 * "RDIdir/RDIprefix*.RDIext". See \ref load_RDIfileList
155 */
156 HCIobservation( const std::string &dir, /**< [in] the directory to search. */
157 const std::string &prefix, /**< [in] the initial part of the file name. Can be empty "".*/
158 const std::string &ext, /**< [in] the extension to append to the file name, must include
159 the '.'.*/
160 const std::string &RDIdir, /**< [in] the directory to search for the reference files.*/
161 const std::string &RDIprefix, /**< [in] the initial part of the file name for the reference files.
162 Can be empty "".*/
163 const std::string &RDIext = "" /**< [in] [optional] the extension to append to the RDI file name,
164 must include the '.'. If empty "" then same
165 extension as target files is used.*/
166 );
167
168 /// Construct using a file containing the target file list and a file containing the RDI target file list
169 /** Populates the \ref m_fileList vector by reading the file, which should be a single
170 * column of new-line delimited file names.
171 *
172 * Populates the \ref m_RDIfileList vector by reading the file, which should be a single
173 * column of new-line delimited file names.
174 */
175 HCIobservation( const std::string &fileListFile, ///< [in] a file name path to read for the target file names.
176 const std::string &RDIfileListFile ///< [in] a file name path to read for the reference file names.
177 );
178
179 /// Load the file list
180 /** Populates the \ref m_fileList vector by searching on disk for files which match the given parameters.
181 * Uses \ref mx::getFileNames to search for all files which match "dir/prefix*.ext".
182 *
183 */
184 void load_fileList( const std::string &dir, ///< [in] the directory to search.
185 const std::string &prefix, ///< [in] the initial part of the file name. Can be empty "".
186 const std::string &ext ///< [in] the extension to append to the file name, which must include
187 ///< the '.'. Can be empty "".
188 );
189
190 /// Load the file list from a file
191 /** Populates the \ref m_fileList vector by reading the file, which should be a single
192 * column of new-line delimited file names.
193 *
194 */
195 void load_fileList( const std::string &fileListFile /**< [in] a file name path to read.*/ );
196
197 /// Load the RDI basis file list
198 /** Populates the \ref m_RDIfileList vector by searching on disk for files which match the given parameters.
199 * Uses \ref mx::getFileNames to search for all files which match "dir/prefix*.ext".
200 *
201 */
202 void load_RDIfileList( const std::string &dir, ///< [in] the directory to search.
203 const std::string &prefix, ///< [in] the initial part of the file name. Can be empty "".
204 const std::string &ext ///< [in] the extension to append to the file name, which must include
205 ///< the '.'. Can be empty "".
206 );
207
208 /// Load the RDI basis file list from a file
209 /** Populates the \ref m_fileList vector by reading the file, which should be a single
210 * column of new-line delimited file names.
211 *
212 */
213 void load_RDIfileList( const std::string &fileListFile /**< [in] a file name path to read.*/ );
214
215 ///@}
216
217 ///\name The Input Target Images
218 /** @{
219 */
220
221 /// The target image cube
223
224 int m_Nims{ 0 }; ///< Number of images in m_tgtIms
225 int m_Nrows{ 0 }; ///< Number of rows of the images in m_tgtIms
226 int m_Ncols{ 0 }; ///< Number of columns of the images in m_tgtIms
227 int m_Npix{ 0 }; ///< Pixels per image, that is Nrows*Ncols
228
229 /// Vector of target image times, in MJD.
230 std::vector<double> m_imageMJD;
231
232 /// Vector of FITS headers, one per file, populated with the values for the keywords.
233 std::vector<fits::fitsHeader> m_heads;
234
235 /// Whether or not the m_fileList has been read.
236 bool m_filesRead{ false };
237
238 /// Whether or not the specified files have been deleted from m_fileList
239 bool m_filesDeleted{ false };
240
241 /// Vector to hold the image weights read from the m_weightFile.
242 /** After readWeights is executed by readFiles, this will contain the normalized weights.
243 * \todo check how comboWeights are handled in coadding
244 */
245 std::vector<realT> m_comboWeights;
246
247 ///@}
248
249 ///\name The Input Reference Images
250 /** @{
251 */
252
253 /// The optional reference image cube
255
256 /// Vector of reference image times, in MJD.
257 std::vector<double> m_RDIimageMJD;
258
259 /// Vector of FITS headers, one per reference file, populated with the values for the keywords.
260 std::vector<fits::fitsHeader> m_RDIheads;
261
262 /// Whether or not the reference files have been read.
263 bool m_RDIfilesRead{ false };
264
265 /// Whether or not the specified files have been deleted from m_RDIfileList
266 bool m_RDIfilesDeleted{ false };
267
268 ///@}
269
270 ///\name The Reduced Data
271 /** @{
272 */
273 /// The PSF subtracted images
274 /** This is a vector of cubes so that it can contain results from different reductions,
275 * e.g. different modes when using KLIP.
276 */
277 std::vector<eigenCube<realT>> m_psfsub;
278
279 /// The final combined images, one for each cube in psfsub.
281
282 ///@}
283
284 /** \name Target File Reading
285 * Options to control which files are read, how they are read, what meta data is extracted
286 * from FITS headers, sizing and masking, etc.
287 * @{
288 */
289
290 /// The list of files to read in.
291 /** This can be set on construction or by calling \ref load_fileList
292 */
293 std::vector<std::string> m_fileList;
294
295 /// Specify how many files from m_fileList to delete from the front of the list
297
298 /// Specify how many files from m_fileList to delete from the back of the list
299 int m_deleteBack{ 0 };
300
301 /// File containing 2 space-delimited columns of fileVame qualityValue pairs.
302 /** If this is not empty and \ref qualityThreshold is > 0, then only images where
303 * qualityValue >= qualityThreshold are read and processed.
304 *
305 * The only restriction on qualityThreshold is that it is > 0. It is intendend to be
306 * something like Strehl ratio.
307 */
308 std::string m_qualityFile;
309
310 /// Threshold to apply to qualityValues read from \ref qualityFile.
311 /** If <= 0, then thresholding is not performed.
312 */
314
315 /// Just prints the names and qualities of the files which pass threshold, and stop.
316 /** Useful mainly for debugging.
317 */
318 bool m_thresholdOnly{ false };
319
320 /// Name of the keyword to use for the image date.
321 /** Specifies the keyword corresponding to the date. This is
322 * the "DATE" keyword for file write time, and usually "DATE-OBS" for actual observation time.
323 *
324 * Default is "DATE-OBS".
325 *
326 * If empty "", then image date is not read.
327 */
328 std::string m_MJDKeyword{ "DATE-OBS" };
329
330 /// Whether or not the date is in ISO 8601 format
331 bool m_MJDisISO8601{ true };
332
333 /// If the date is not ISO 8601, this specifies the conversion to Julian Days (i.e. seconds to days)
335
336 /// Vector of FITS header keywords to read from the files in m_fileList.
337 std::vector<std::string> m_keywords;
338
339 /// Set the image size. Images are cut down to this size after reading.
340 /** Set to <= 0 to use images uncut.
341 *
342 * Image sizes are not increased if this is larger than their size on disk.
343 */
344 int m_imSize{ 0 };
345
346 /// Read the list of files, cut to size, and preprocess.
347 /**
348 * \returns 0 on success, -1 on error.
349 */
350 int readFiles();
351
352 /// Perform post-read actions for the target images, for use by derived classes
353 virtual int postReadFiles();
354
355 /// Perform post-coadd actions for the target images, for use by derived classes.
356 /**
357 * \returns 0 on success
358 * \returns <0 on error.
359 */
360 virtual int postCoadd();
361
362 ///@}
363
364 /** \name Reference File Reading
365 * For RDI, Options to control which files are read, how they are read, what meta data is extracted
366 * from FITS headers, sizing and masking, etc.
367 * @{
368 */
369
370 /// The list of files to read in for the RDI basis.
371 /** This is set by calling \ref load_RDIfileList
372 */
373 std::vector<std::string> m_RDIfileList;
374
375 /// Specify how many files from m_RDIfileList to delete from the front of the list
377
378 /// Specify how many files from m_RDIfileList to delete from the back of the list
380
381 /// File containing 2 space-delimited columns of fileMame qualityValue pairs for the reference images.
382 /** If this is not empty and \ref m_RDIqualityThreshold is > 0, then only images where
383 * qualityValue >= qualityThreshold are read and processed.
384 *
385 * The only restriction on m_RDIqualityThreshold is that it is > 0. It is intendend to be
386 * something like Strehl ratio.
387 */
388 std::string m_RDIqualityFile;
389
390 /// Threshold to apply to qualityValues read from \ref qualityFile.
391 /** If <= 0, then thresholding is not performed.
392 */
394
395 /// Vector of FITS header keywords to read from the files in m_fileList.
396 std::vector<std::string> m_RDIkeywords;
397
398 /// Read the list of reference files, cut to size, and preprocess.
399 /** The target files must be read with \ref readFiles() before calling this method.
400 *
401 * \returns 0 on success, -1 on error.
402 */
403 int readRDIFiles();
404
405 /// Perform post-read actions for the RDI images, for use by derived classes
406 virtual int postRDIReadFiles();
407
408 /// Perform post-coadd actions, for use by derived classes.
409 /** A key example is to update keywords after the averaging occurs in coaddImages().
410 *
411 * \returns 0 on success
412 * \returns <0 on error.
413 */
414 virtual int postRDICoadd();
415
416 ///@}
417 //--RDI
418
419 /** \name Thresholding
420 * Thresholds are applied to a list of files before it is read, based on the image qualities supplied.
421 * @{
422 */
423
424 /// Read the image qualities from a qualityFile and apply the threshold to a fileList
425 /** This is called by readFiles().
426 *
427 * \returns 0 on success, -1 on error.
428 */
429 int threshold( std::vector<std::string> &fileList, ///< [in.out] the fileList to threshold
430 const std::string &qualityFile, ///< [in] the path to the file containing qualities, one per file.
431 realT qualityThreshold ///< [in] the quality threshold to apply
432 );
433
434 ///@}
435
436 /** \name Coadding
437 * These parameters control whether and how the images are coadded after being read. Coadding can
438 * be done up to a given number of images, and/or a given elapsed time.
439 *
440 * Averages the values of given Keywords as well.
441 * @{
442 */
443
444 /// Determine how to coadd the raw images.
445 /** Possibilities are
446 * - HCI::noCombine -- [default] do not combine. This turns off coadding.
447 * - HCI::medianCombine -- coadded image is the median
448 * - HCI::meanCombine -- coadded image is the simple mean
449 *
450 * No other types of combination are currently supported for coadding.
451 */
453
454 /// Maximum number of images to coadd
456
457 /// Maximum elapsed time over which to coadd the images.
459
460 /// The values of these keywords will be averaged and replaced.
461 std::vector<std::string> m_coaddKeywords;
462
463 /// Coadd the images
464 void coaddImages( int coaddCombineMethod,
465 int coaddMaxImno,
466 int coaddMaxTime,
467 std::vector<std::string> &coaddKeywords,
468 std::vector<double> &imageMJD,
469 std::vector<fits::fitsHeader> &heads,
470 eigenCube<realT> &ims );
471
472 ///@} -- coadding
473
474 /** \name Masking
475 * A 1/0 mask can be supplied, which is used in pre-processing and in image combination.
476 * @{
477 */
478
479 /// Specify a mask file to apply
480 /**No mask is applied if this is empty.
481 */
482 std::string m_maskFile;
483
484 eigenImageT m_mask; ///< The mask
485
486 eigenCube<realT> m_maskCube; /**< A cube of masks, one for each input image, which may be modified
487 versions (e.g. rotated) of mask. */
488
489 /// Read the mask file, resizing to imSize if needed.
490 void readMask();
491
492 /// Populate the mask cube which is used for post-processing.
493 /** Derived classes can do this as appropriate, e.g. by rotating the mask.
494 * \throws mx::err::invalidconfig if mask is not the same size as the images
495 */
496 virtual void makeMaskCube();
497
498 ///@}
499
500 public:
501 /** \name Pre-Processing
502 * These options control the pre-processing masking and filtering.
503 * They are performed in the following order:
504 * -# mask applied (enabled by m_preProcess_mask
505 * -# radial profile subtraction (enabled by m_preProcess_subradprof)
506 * -# mask applied (enabled by m_preProcess_mask
507 * -# symmetric median unsharp mask (m_preProcess_gaussUSM_fwhm)
508 * -# symmetric Gaussian unsharp mask (m_preProcess_gaussUSM_fwhm)
509 * -# mask applied (enabled by m_preProcess_mask
510 * -# azimuthal unsharp mask (m_preProcess_azUSM_azW, and m_preProcess_azUSM_radW)
511 * -# mask applied (enabled by m_preProcess_mask)
512 * -# mean subtraction (enabled by m_preProcess_meanSubMethod)
513 * -# mask applied (enabled by m_preProcess_mask)
514 * -# pixel time-series normalization (enabled by m_preProcess_pixelTSNormMethod)
515 * @{
516 */
517
518 bool m_skipPreProcess{ false }; ///< Don't do any of the pre-processing steps (including coadding).
519
520 bool m_preProcess_beforeCoadd{ false }; ///< controls whether pre-processing takes place before or after coadding
521
522 bool m_preProcess_mask{ true }; ///< If true, the mask is applied during each pre-processing step.
523
524 bool m_preProcess_subradprof{ false }; ///< If true, a radial profile is subtracted from each image.
525
526 /// Azimuthal boxcar width for azimuthal unsharp mask [pixels]
527 /** If this is 0 then azimuthal-USM is not performed.
528 */
530
531 /// Mazimum azimuthal boxcar width for azimuthal unsharp mask [degrees]
532 /** Limits width close to center, preventing wrap-around. Default is 45 degrees. Set to 0 for no maximum.
533 */
535
536 /// Radial boxcar width for azimuthal unsharp mask [pixels]
537 /** If this is 0 then azimuthal-USM is not performed.
538 */
540
541 /// Kernel FWHM for symmetric box median unsharp mask (USM)
542 /** USM is not performed if this is 0.
543 */
545
546 /// Kernel FWHM for symmetric Gaussian unsharp mask (USM)
547 /** USM is not performed if this is 0.
548 */
550
551 /// The mean subtraction method during pre-processing
552 /** Can only be none, meanImage, or meadianImage
553 */
555
556 /// Specify if each pixel time-series is normalized
557 /** This normalizaton is applied after centering. Can have the following values:
558 * - <b>HCI::pixelTSNormMethod::none</b>: no normalization (the default)
559 * - <b>HCI::pixelTSNormMethod::rms</b>: divide by the time-series rms
560 * - <b>HCI::pixelTSNormMethod::rmsSigmaClipped</b>: divide by the sigma-slipped time-series rms.
561 * The sigma is provided by m_preProcess_pixelTSSigma.
562 */
564
565 realT m_pixelTSSigma{ 3 }; ///< Sigma-clipping parameter for pixel time-series normalization
566
567 /// Set path and file prefix to output the pre-processed images.
568 /** If empty, then pre-processed images are not output.
569 */
571
572 /// If true, then we stop after pre-processing.
573 bool m_preProcess_only{ false };
574
575 /// Do the pre-processing
576 void preProcess( eigenCube<realT> &ims /**< [in] the image cube, should be either m_tgtIms or m_refIms */ );
577
578 /// Do mean subtraction as part of pre-processing
579 void preProcess_meanSub( eigenCube<realT> &ims /**< [in] the image cube, should be either m_tgtIms or m_refIms */ );
580
581 /// Do pixel time-series normalization as part of pre-processing
582 void
583 preProcess_pixelTSNorm( eigenCube<realT> &ims /**< [in] the image cube, should be either m_tgtIms or m_refIms */ );
584
585 ///@}
586
587 /** \name Image Combination
588 * These options control how the final image combination is performed.
589 * @{
590 */
591
592 /// Determine how to combine the PSF subtracted images
593 /** Possibilities are
594 * - HCI::noCombine -- do not combine
595 * - HCI::medianCombine -- [default] final image is the median
596 * - HCI::meanCombine -- final image is the simple mean
597 * - HCI::weightedMeanCombine -- final image is the weighted mean. m_weightFile must be provided.
598 * - HCI::sigmaMeanCombine -- final image is sigma clipped mean. If m_sigmaThreshold <= 0, then it reverts to
599 * meanCombine.
600 */
602
603 /// Specifies a file containing the image weights, for combining with weighted mean.
604 /** This 2-column space-delimited ASCII file containing filenames and weights. It must be specified before
605 * readFiles() is executed. In readFiles this is loaded after the m_fileList is cutdown and matched to the
606 * remaining files.
607 */
608 std::string m_weightFile;
609
610 /// The standard deviation threshold used if combineMethod == HCI::sigmaMeanCombine.
612
613 /// The minimum fraction of good (un-masked) pixels to include in the final combination (0.0 to 1.0). If not met,
614 /// then the pixel will be NaN-ed.
616
617 /// Read the image weights from m_weightFile
618 /** This is called by readFiles().
619 *
620 * \returns 0 on success, -1 on error.
621 */
622 int readWeights();
623
624 /// Combine the images into a single final image.
625 /** Images are combined by the method specified in \ref combineMethod
626 */
627 void combineFinim();
628
629 ///@}
630
631 /** \name Output
632 * These options control the ouput of the final combined images and the individual PSF subtracted images.
633 * @{
634 */
635
636 /// Location for temporary auxilliary output files (e.g. masks)
637 std::string m_auxDataDir{ "/tmp/klipReduceAux/" };
638
639 /// Whether or not to move the temp. aux files.
640 bool m_moveAuxDataDir{ true };
641
642 /// Set whether the final combined image is written to disk
644
645 /// The directory where to write output files.
646 std::string m_outputDir;
647
648 /// The base file name of the output final image
649 /** The complete name is formed by combining with a sequential number and the ".fits" extension.
650 * that is: m_finimName0000.fits. This behavior can be modified with m_exactFinimName.
651 */
652 std::string m_finimName{ "finim_" };
653
654 /// Use m_finimName exactly as specified, without appending a number or an extension.
655 /** Output is still FITS format, regardless of extension. This will overwrite
656 * an existing file without asking.
657 */
658 bool m_exactFinimName{ false };
659
660 /// Controls whether or not the individual PSF subtracted images are written to disk.
661 /** - true -- write to disk
662 * - false -- [default] don't write to disk
663 */
664 bool m_doOutputPSFSub{ false };
665
666 /// Prefix of the FITS file names used to write individual PSF subtracted images to disk if m_doOutputPSFSub is
667 /// true.
668 std::string m_PSFSubPrefix;
669
670 /// Output the pre-processed target images
671 void outputPreProcessed();
672
673 /// Output the pre-processed reference images
675
676 /// Fill in the HCIobservation standard FITS header
677 /**
678 */
679 void stdFitsHeader( fits::fitsHeader &head /**< [in.out] the fistHeader structure which will
680 have cards appended to it. */
681 );
682
683 /// Write the final combined image to disk
684 /**
685 */
686 void writeFinim( fits::fitsHeader *addHead = 0 );
687
688 /// Write the PSF subtracted images to disk
689 /**
690 */
691 void outputPSFSub( fits::fitsHeader *addHead = 0 );
692
693 ///@}
694
695 /// Read in already PSF-subtracted files
696 /** Used to take up final processing after applying some non-klipReduce processing steps to
697 * PSF-subtracted images.
698 */
699 int readPSFSub( const std::string &dir, const std::string &prefix, const std::string &ext, size_t nReductions );
700
701 double t_begin{ 0 };
702 double t_end{ 0 };
703
704 double t_load_begin{ 0 };
705 double t_load_end{ 0 };
706
707 double t_coadd_begin{ 0 };
708 double t_coadd_end{ 0 };
709
710 double t_preproc_begin{ 0 };
711 double t_preproc_end{ 0 };
712
713 double t_azusm_begin{ 0 };
714 double t_azusm_end{ 0 };
715
716 double t_gaussusm_begin{ 0 };
717 double t_gaussusm_end{ 0 };
718
719 double t_combo_begin{ 0 };
720 double t_combo_end{ 0 };
721};
722
723// -- construction and initialization
724
725template <typename _realT>
729
730template <typename _realT>
731HCIobservation<_realT>::HCIobservation( const std::string &dir, const std::string &prefix, const std::string &ext )
732{
733 load_fileList( dir, prefix, ext );
734}
735
736template <typename _realT>
737HCIobservation<_realT>::HCIobservation( const std::string &fileListFile )
738{
739 load_fileList( fileListFile );
740}
741
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 )
749{
750 std::cerr << "HCI 6\n";
751
752 load_fileList( dir, prefix, ext );
753
754 std::string re;
755 if( RDIext == "" )
756 re = ext;
757 else
758 re = RDIext;
759
760 load_RDIfileList( RDIdir, RDIprefix, re );
761}
762
763template <typename _realT>
764HCIobservation<_realT>::HCIobservation( const std::string &fileListFile, const std::string &RDIfileListFile )
765{
766 load_fileList( fileListFile );
767 load_RDIfileList( RDIfileListFile );
768}
769
770template <typename _realT>
771void HCIobservation<_realT>::load_fileList( const std::string &dir, const std::string &prefix, const std::string &ext )
772{
773 m_fileList = ioutils::getFileNames( dir, prefix, "", ext );
774
775 m_filesDeleted = false;
776}
777
778template <typename _realT>
779void HCIobservation<_realT>::load_fileList( const std::string &fileListFile )
780{
781 ioutils::readColumns( fileListFile, m_fileList );
782 m_filesDeleted = false;
783}
784
785template <typename _realT>
786void HCIobservation<_realT>::load_RDIfileList( const std::string &dir,
787 const std::string &prefix,
788 const std::string &ext )
789{
790 m_RDIfileList = ioutils::getFileNames( dir, prefix, "", ext );
791 m_RDIfilesDeleted = false;
792}
793
794template <typename _realT>
795void HCIobservation<_realT>::load_RDIfileList( const std::string &fileListFile )
796{
797 ioutils::readColumns( fileListFile, m_RDIfileList );
798 m_RDIfilesDeleted = false;
799}
800
801// --< construction and initialization
802
803template <typename realT>
805{
806 if( m_fileList.size() == 0 )
807 {
808 mxThrowException( err::invalidconfig,
809 "HCIobservation<realT>::readFiles",
810 "The target fileList has 0 length, there are no files to be read." );
811 }
812
813 // First make the list deletions
814 if( !m_filesDeleted )
815 {
816 if( m_deleteFront > 0 )
817 {
818 m_fileList.erase( m_fileList.begin(), m_fileList.begin() + m_deleteFront );
819 }
820
821 if( m_deleteBack > 0 )
822 {
823 m_fileList.erase( m_fileList.end() - m_deleteBack, m_fileList.end() );
824 }
825 m_filesDeleted = true;
826 }
827
828 if( m_fileList.size() == 0 )
829 {
830 mxThrowException( err::invalidconfig,
831 "HCIobservation<realT>::readFiles",
832 "The target fileList has 0 length, there are no files to be read after deletions." );
833 }
834
835 if( m_qualityFile != "" )
836 {
837 std::cerr << "Thresholding target images...";
838 size_t origsize = m_fileList.size();
839
840 if( threshold( m_fileList, m_qualityFile, m_qualityThreshold ) < 0 )
841 return -1;
842
843 if( m_fileList.size() == 0 )
844 {
845 mxThrowException( err::invalidconfig,
846 "HCIobservation<realT>::readFiles",
847 "The fileList has 0 length, there are no files to be read after thresholding." );
848 }
849
850 std::cerr << "Done. Selected " << m_fileList.size() << " out of " << origsize << "\n";
851
852 if( m_thresholdOnly )
853 {
854 std::cout << "#Files which passed thresholding:\n";
855 for( size_t i = 0; i < m_fileList.size(); ++i )
856 {
857 std::cout << m_fileList[i] << "\n";
858 }
859
860 exit( 0 );
861 }
862 }
863
864 if( m_weightFile != "" )
865 {
866 if( readWeights() < 0 )
867 return -1;
868 }
869
870 /*----- Append the HCI keywords to propagate them if needed -----*/
871
872 fits::fitsHeader head;
873
874 if( m_MJDKeyword != "" )
875 head.append( m_MJDKeyword );
876
877 for( size_t i = 0; i < m_keywords.size(); ++i )
878 {
879 if( head.count( m_keywords[i] ) == 0 )
880 head.append( m_keywords[i] );
881 }
882
883 m_heads.clear(); // This is necessary to make sure heads.resize() copies head on a 2nd call
884 m_heads.resize( m_fileList.size(), head );
885
886 /*----- Read in first image to test size -----*/
887
888 Eigen::Array<realT, Eigen::Dynamic, Eigen::Dynamic> im;
889
890 fits::fitsFile<realT> f( m_fileList[0] );
891
892 f.read( im );
893
894 // We set imSize to match the first image, but we make it a square.
895 if( m_imSize == 0 )
896 {
897 m_imSize = im.rows();
898 if( m_imSize > im.cols() )
899 {
900 m_imSize = im.cols();
901 }
902 }
903 else
904 {
905 // Now make sure we don't read too much.
906 if( m_imSize > im.rows() )
907 {
908 m_imSize = im.rows();
909 }
910 if( m_imSize > im.cols() )
911 {
912 m_imSize = im.cols();
913 }
914 }
915
916 // And now set the read size so we only read what we want.
917 // the +0.1 is just to make sure we don't have a problem with precision (we shouldn't)/
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 ),
920 m_imSize,
921 m_imSize );
922
923 im.resize( m_imSize, m_imSize );
924
925 std::cerr << "image size is " << m_imSize << "\n";
926
927 /**** Right here is we go to coadd.
928 */
929
930 // But if not we just do it
931 m_tgtIms.resize( im.rows(), im.cols(), m_fileList.size() );
932
933 t_load_begin = sys::get_curr_time();
934
935 f.read( m_tgtIms.data(), m_heads, m_fileList );
936
937 f.setReadSize();
938
939 /* read in the image timestamps, depending on how MJD is stored in the header */
940 if( m_MJDKeyword != "" )
941 {
942 m_imageMJD.resize( m_heads.size() );
943
944 if( m_MJDisISO8601 )
945 {
946 for( size_t i = 0; i < m_imageMJD.size(); ++i )
947 {
948 m_imageMJD[i] = sys::ISO8601date2mjd( m_heads[i][m_MJDKeyword].String() );
949 }
950 }
951 else
952 {
953 for( size_t i = 0; i < m_imageMJD.size(); ++i )
954 {
955 m_imageMJD[i] = m_heads[i][m_MJDKeyword].template value<realT>() * m_MJDUnits;
956 }
957 }
958 }
959
960 t_load_end = sys::get_curr_time();
961
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();
966
967 std::cerr << "loading complete\n";
968
969 std::cerr << "zero-ing NaNs\n";
970 zeroNaNCube( m_tgtIms );
971 std::cerr << "done\n";
972
973 /*** Now do the post-read actions ***/
974 if( postReadFiles() < 0 )
975 {
976 return -1;
977 }
978
979 /*** Read in the mask if present ***/
980 readMask();
981
982 /*** Now begin processing ***/
983 if( !m_skipPreProcess )
984 {
985 /*** Now do any pre-processing ***/
986 if( m_preProcess_beforeCoadd )
987 {
988 preProcess( m_tgtIms );
989 }
990
991 if( m_coaddCombineMethod != HCI::noCombine )
992 {
993 std::cerr << "Coadding target images...\n";
994 coaddImages( m_coaddCombineMethod,
995 m_coaddMaxImno,
996 m_coaddMaxTime,
997 m_coaddKeywords,
998 m_imageMJD,
999 m_heads,
1000 m_tgtIms );
1001
1002
1003
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();
1008
1009 std::cerr << "number of target images after coadding: " << m_Nims << "\n";
1010
1011 if( postCoadd() < 0 )
1012 {
1013 std::cerr << "Post coadd error " << __FILE__ << " " << __LINE__ << "\n";
1014 return -1;
1015 }
1016 std::cerr << "Done.\n";
1017
1018 // Re-make the mask cube if we coadded...
1019 if( m_maskFile != "" )
1020 {
1021 makeMaskCube();
1022 }
1023 }
1024
1025 /*** Now do any pre-processing if not done already***/
1026 if( !m_preProcess_beforeCoadd )
1027 {
1028 preProcess( m_tgtIms );
1029 }
1030
1031 outputPreProcessed();
1032 }
1033
1034 m_filesRead = true;
1035
1036 return 0;
1037} // readFiles()
1038
1039template <typename _realT>
1041{
1042 return 0;
1043}
1044
1045template <typename _realT>
1047{
1048 return 0;
1049}
1050
1051//------------------- readRDIFiles
1052template <typename _realT>
1054{
1055
1056 /* First check if the target files have been read */
1057 if( m_Nrows == 0 || m_Ncols == 0 )
1058 {
1059 mxError( "HCIobservation", MXE_PARAMNOTSET, "The target image size must be set before reading RDI files." );
1060 return -1;
1061 }
1062
1063 if( m_RDIfileList.size() == 0 )
1064 {
1065 mxError( "HCIobservation", MXE_FILENOTFOUND, "The RDI fileList has 0 length, there are no files to be read." );
1066 return -1;
1067 }
1068
1069 // First make the list deletions
1070 if( !m_RDIfilesDeleted )
1071 {
1072 if( m_RDIdeleteFront > 0 )
1073 {
1074 m_RDIfileList.erase( m_RDIfileList.begin(), m_RDIfileList.begin() + m_RDIdeleteFront );
1075 }
1076
1077 if( m_RDIdeleteBack > 0 )
1078 {
1079 m_RDIfileList.erase( m_RDIfileList.end() - m_RDIdeleteBack, m_RDIfileList.end() );
1080 }
1081 m_RDIfilesDeleted = true;
1082 }
1083
1084 if( m_RDIfileList.size() == 0 )
1085 {
1086 mxError( "HCIobservation",
1087 MXE_FILENOTFOUND,
1088 "The RDI fileList has 0 length, there are no files to be read after deletions." );
1089 return -1;
1090 }
1091
1092 if( m_RDIqualityFile != "" )
1093 {
1094 std::cerr << "Thresholding RDI images...";
1095 size_t origsize = m_RDIfileList.size();
1096
1097 if( threshold( m_RDIfileList, m_RDIqualityFile, m_RDIqualityThreshold ) < 0 )
1098 return -1;
1099
1100 if( m_RDIfileList.size() == 0 )
1101 {
1102 mxError( "HCIobservation",
1103 MXE_FILENOTFOUND,
1104 "The fileList has 0 length, there are no files to be read after thresholding." );
1105 return -1;
1106 }
1107
1108 std::cerr << "Done. Selected " << m_RDIfileList.size() << " out of " << origsize << "\n";
1109 }
1110
1111 /*----- Append the HCI keywords to propagate them if needed -----*/
1112
1113 fits::fitsHeader head;
1114
1115 if( m_MJDKeyword != "" )
1116 head.append( m_MJDKeyword ); // Currently assuming the MJD keyword will be the same
1117
1118 for( size_t i = 0; i < m_RDIkeywords.size(); ++i )
1119 {
1120 head.append( m_RDIkeywords[i] );
1121 }
1122
1123 m_RDIheads.clear(); // This is necessary to make sure heads.resize() copies head on a 2nd call
1124 m_RDIheads.resize( m_RDIfileList.size(), head );
1125
1126 /*----- Read in first image to adjust size ----*/
1127 Eigen::Array<realT, Eigen::Dynamic, Eigen::Dynamic> im;
1128
1129 fits::fitsFile<realT> f( m_RDIfileList[0] );
1130
1131 f.read( im );
1132
1133 if( im.rows() < m_imSize || im.cols() < m_imSize )
1134 {
1135 mxError( "HCIobservation", MXE_SIZEERR, "The reference images are too small, do not match the target images." );
1136 return -1;
1137 }
1138
1139 // And now set the read size so we only read what we want.
1140 // the +0.1 is just to make sure we don't have a problem with precision (we shouldn't)/
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 ),
1143 m_imSize,
1144 m_imSize );
1145
1146 m_refIms.resize( m_imSize, m_imSize, m_RDIfileList.size() );
1147
1148 t_load_begin = sys::get_curr_time();
1149
1150 f.read( m_refIms.data(), m_RDIheads, m_RDIfileList );
1151
1152 f.setReadSize();
1153
1154 if( m_MJDKeyword != "" )
1155 {
1156 m_RDIimageMJD.resize( m_RDIheads.size() );
1157
1158 if( m_MJDisISO8601 )
1159 {
1160 for( size_t i = 0; i < m_RDIimageMJD.size(); ++i )
1161 {
1162 m_RDIimageMJD[i] = sys::ISO8601date2mjd( m_RDIheads[i][m_MJDKeyword].String() );
1163 }
1164 }
1165 else
1166 {
1167 for( size_t i = 0; i < m_RDIimageMJD.size(); ++i )
1168 {
1169 m_RDIimageMJD[i] = m_RDIheads[i][m_MJDKeyword].template value<realT>() * m_MJDUnits;
1170 }
1171 }
1172 }
1173
1174 t_load_end = sys::get_curr_time();
1175
1176 std::cerr << "loading complete\n";
1177
1178 std::cerr << "zero-ing NaNs\n";
1179 zeroNaNCube( m_refIms );
1180 std::cerr << "Done.\n";
1181
1182 /*** Now do the post-read actions ***/
1183 if( postRDIReadFiles() < 0 )
1184 return -1;
1185
1186 /*** Now begin processing ***/
1187 if( !m_skipPreProcess )
1188 {
1189 /*** Now do any pre-processing ***/
1190 if( m_preProcess_beforeCoadd )
1191 {
1192 preProcess( m_refIms );
1193 }
1194
1195 if( m_coaddCombineMethod != HCI::noCombine )
1196 {
1197 std::cerr << "Coadding reference images...\n";
1198 coaddImages( m_coaddCombineMethod,
1199 m_coaddMaxImno,
1200 m_coaddMaxTime,
1201 m_coaddKeywords,
1202 m_RDIimageMJD,
1203 m_RDIheads,
1204 m_refIms );
1205
1206 std::cerr << "number of reference images after coadding: " << m_refIms.planes() << "\n";
1207
1208 if( postRDICoadd() < 0 )
1209 {
1210 std::cerr << "Post coadd error " << __FILE__ << " " << __LINE__ << "\n";
1211 return -1;
1212 }
1213 std::cerr << "Done.\n";
1214 }
1215
1216 /*** Now do any pre-processing if not done already***/
1217 if( !m_preProcess_beforeCoadd )
1218 {
1219 preProcess( m_refIms );
1220 }
1221
1222 // outputRDIPreProcessed();
1223 }
1224
1225 m_RDIfilesRead = true;
1226
1227 return 0;
1228} // readRDIFiles()
1229
1230template <typename _realT>
1232{
1233 return 0;
1234}
1235
1236template <typename _realT>
1238{
1239 return 0;
1240}
1241
1242template <typename _realT>
1243int HCIobservation<_realT>::threshold( std::vector<std::string> &fileList,
1244 const std::string &qualityFile,
1245 realT qualityThreshold )
1246{
1247 if( qualityFile == "" )
1248 {
1249 mxError( "HCIobservation::threshold", MXE_PARAMNOTSET, "qualityFile not set" );
1250 return -1;
1251 }
1252
1253 int origsize = fileList.size();
1254
1255 std::vector<std::string> qfileNames;
1256 std::vector<realT> imQ;
1257
1258 // Read the quality file and load it into a map
1259 ioutils::readColumns( qualityFile, qfileNames, imQ );
1260
1261 std::map<std::string, realT> quality;
1262 for( size_t i = 0; i < qfileNames.size(); ++i )
1263 quality[ioutils::pathFilename( qfileNames[i].c_str() )] = imQ[i];
1264
1265 realT q;
1266
1267 for( size_t i = 0; i < fileList.size(); ++i )
1268 {
1269 try
1270 {
1271 q = quality.at( ioutils::pathFilename( fileList[i].c_str() ) );
1272 }
1273 catch( ... )
1274 {
1275 q = qualityThreshold - 1; // Cause it to be erased
1276 }
1277
1278 if( q < qualityThreshold )
1279 {
1280 fileList.erase( fileList.begin() + i );
1281 --i;
1282 }
1283 }
1284
1285 return 0;
1286}
1287
1288template <typename _realT>
1289void HCIobservation<_realT>::coaddImages( int coaddCombineMethod,
1290 int coaddMaxImno,
1291 int coaddMaxTime,
1292 std::vector<std::string> &coaddKeywords,
1293 std::vector<double> &imageMJD,
1294 std::vector<fits::fitsHeader> &heads,
1295 eigenCube<realT> &ims )
1296{
1297 std::cerr << "***************************************************************\n";
1298 std::cerr << " *** WARNING *** \n";
1299 std::cerr << " coadding is poorly tested. proceed with caution. \n";
1300 std::cerr << "***************************************************************\n";
1301
1302 // Validate setup
1303 if( coaddMaxImno <= 0 && coaddMaxTime <= 0 )
1304 return;
1305
1306 // Validate combine method
1307 if( coaddCombineMethod == HCI::noCombine )
1308 return;
1309
1310 int Nims = ims.planes();
1311 int Nrows = ims.rows();
1312 int Ncols = ims.cols();
1313
1314 t_coadd_begin = sys::get_curr_time();
1315
1316 std::vector<eigenImageT> coadds;
1317 std::vector<std::vector<std::string>> coaddFileNames;
1318
1319 // We do all math here in double, to avoid losing precision
1320 std::vector<double> avgMJD;
1321 std::vector<double> startMJD;
1322 std::vector<double> endMJD;
1323
1324 std::vector<std::vector<double>> avgVals;
1325 std::vector<std::vector<double>> startVals;
1326 std::vector<std::vector<double>> endVals;
1327
1328 int combineMethod = HCI::medianCombine;
1329 if( coaddCombineMethod == HCI::meanCombine )
1330 combineMethod = HCI::meanCombine;
1331
1332 // Index range of images for next coadd
1333 int im0, imF;
1334 im0 = 0;
1335 imF = 1;
1336
1337 // Accumulate images to coadd into a cube
1338 eigenCube<realT> imsToCoadd;
1339
1340 // The filenames of the images included in the coadd
1341 std::vector<std::string> imsCoadded;
1342
1343 // Temporary for combination.
1344 eigenImageT coadd;
1345
1346 // Accumulate values
1347 double initMJD;
1348
1349 std::vector<double> initVals;
1350 initVals.resize( coaddKeywords.size() );
1351
1352 std::vector<double> startVal;
1353 startVal.resize( coaddKeywords.size() );
1354
1355 std::vector<double> endVal;
1356 endVal.resize( coaddKeywords.size() );
1357
1358 while( im0 < Nims )
1359 {
1360 // Initialize the accumulators
1361 initMJD = imageMJD[im0];
1362 startMJD.push_back( initMJD );
1363 endMJD.push_back( initMJD );
1364
1365 for( size_t i = 0; i < coaddKeywords.size(); ++i )
1366 {
1367 initVals[i] = heads[im0][coaddKeywords[i]].value<double>();
1368 startVal[i] = initVals[i];
1369 endVal[i] = initVals[i];
1370 }
1371
1372 // Now increment imF, then test whether each variable is now outside the range
1373 bool increment = true;
1374 while( increment == true )
1375 {
1376 ++imF;
1377
1378 if( imF >= Nims ) // out of images
1379 {
1380 imF = Nims;
1381 increment = false;
1382 break;
1383 }
1384
1385 if( imF - im0 > coaddMaxImno && coaddMaxImno > 0 ) // too many images in coadd
1386 {
1387 --imF;
1388 increment = false;
1389 break;
1390 }
1391
1392 ///\todo this should really include end of exposure too.
1393 if( ( imageMJD[imF] - imageMJD[im0] ) * 86400.0 > coaddMaxTime && coaddMaxTime > 0 )
1394 {
1395 --imF;
1396 increment = false;
1397 break;
1398 }
1399
1400 } // while(increment == true)
1401 // At this point, imF is the first image NOT to include in the next coadd.
1402
1403 // Now actually do the accumulation
1404 ///\todo this needs to handle averaging of angles
1405 for( int imno = im0 + 1; imno < imF; ++imno )
1406 {
1407 initMJD += imageMJD[imno];
1408 endMJD.back() = imageMJD[imno]; // After the last one, this will be the last one
1409
1410 for( size_t i = 0; i < coaddKeywords.size(); ++i )
1411 {
1412 endVal[i] =
1413 heads[imno][coaddKeywords[i]].value<double>(); // After the last one, this will be the last one
1414 initVals[i] += endVal[i];
1415 }
1416 }
1417
1418 // And then turn them into an average
1419 initMJD /= ( imF - im0 );
1420 for( size_t i = 0; i < coaddKeywords.size(); ++i )
1421 {
1422 initVals[i] /= ( imF - im0 );
1423 }
1424
1425 // Extract the images into the temporary
1426 imsToCoadd.resize( Nrows, Ncols, imF - im0 );
1427 imsCoadded.resize( imF - im0 );
1428 for( int i = 0; i < ( imF - im0 ); ++i )
1429 {
1430 imsToCoadd.image( i ) = ims.image( im0 + i );
1431 imsCoadded[i] = ioutils::pathFilename( m_fileList[im0 + i] );
1432 }
1433
1434 // Here do the combine and insert into the vector
1435 if( combineMethod == HCI::medianCombine )
1436 {
1437 imsToCoadd.median( coadd );
1438 }
1439 if( combineMethod == HCI::meanCombine )
1440 {
1441 imsToCoadd.mean( coadd );
1442 }
1443 coadds.push_back( coadd );
1444 coaddFileNames.push_back( imsCoadded );
1445
1446 // Insert the new averages
1447 avgMJD.push_back( initMJD );
1448 avgVals.push_back( initVals );
1449 startVals.push_back( startVal );
1450 endVals.push_back( endVal );
1451
1452 im0 = imF;
1453 imF = im0 + 1;
1454 } // while(im0 < Nims)
1455
1456 // Now resize ims and copy the coadds to the new cube
1457 ims.resize( Nrows, Ncols, coadds.size() );
1458 Nims = coadds.size();
1459
1460 for( int i = 0; i < Nims; ++i )
1461 {
1462 ims.image( i ) = coadds[i];
1463 }
1464
1465 // Now deal with imageMJD and headers
1466 imageMJD.erase( imageMJD.begin() + Nims, imageMJD.end() );
1467 heads.erase( heads.begin() + Nims, heads.end() );
1468
1469 for( int i = 0; i < Nims; ++i )
1470 {
1471 imageMJD[i] = avgMJD[i];
1472 heads[i][m_MJDKeyword].value( mx::sys::ISO8601DateTimeStrMJD( imageMJD[i], 1 ) );
1473 heads[i]["START " + m_MJDKeyword].value( mx::sys::ISO8601DateTimeStrMJD( startMJD[i], 1 ) );
1474 heads[i]["END " + m_MJDKeyword].value( mx::sys::ISO8601DateTimeStrMJD( endMJD[i], 1 ) );
1475 heads[i].append( "DELTA " + m_MJDKeyword,
1476 ( endMJD[i] - startMJD[i] ) * 86400,
1477 "change in " + m_MJDKeyword + " in seconds." );
1478
1479 for( size_t j = 0; j < coaddKeywords.size(); ++j )
1480 {
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] );
1485 }
1486
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 )
1489 {
1490 // fits::fitsHistoryType t;
1491 // fits::fitsHeaderCard c("HISTORY", t , "coadded file: " + coaddFileNames[i][j]);
1492 heads[i].append( "HISTORY", fits::fitsHistoryType(), "coadded file: " + coaddFileNames[i][j] );
1493 }
1494 }
1495
1496 t_coadd_end = sys::get_curr_time();
1497
1498} // void HCIobservation<_realT>::coaddImages()
1499
1500template <typename _realT>
1502{
1503 /*** Load the mask ***/
1504 if( m_maskFile != "" )
1505 {
1506 std::cerr << "creating mask cube\n";
1508 ff.read( m_mask, m_maskFile );
1509
1510 ///\todo here re-size mask if needed to match imSize
1511 if( m_mask.rows() > m_imSize || m_mask.cols() > m_imSize )
1512 {
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 ) ),
1515 m_imSize,
1516 m_imSize );
1517 m_mask = tmask;
1518 }
1519
1520 makeMaskCube();
1521 }
1522}
1523
1524template <typename realT>
1526{
1527 if( m_mask.rows() != m_Nrows || m_mask.cols() != m_Ncols )
1528 {
1529 // clang-format off
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";
1535 // clang-format on
1536
1537 mxThrowException( err::invalidconfig, "HCIobservation<realT>::makeMaskCube", message );
1538 }
1539
1540 m_maskCube.resize( m_Nrows, m_Ncols, m_Nims );
1541
1542 for( int i = 0; i < m_Nims; ++i )
1543 {
1544 m_maskCube.image( i ) = m_mask;
1545 }
1546}
1547
1548template <typename _realT>
1550{
1551 t_preproc_begin = sys::get_curr_time();
1552
1553 // The mask is applied first, and then after each subsequent P.P. step.
1554 if( m_maskFile != "" && m_preProcess_mask )
1555 {
1556 std::cerr << "Masking . . .\n";
1557#pragma omp parallel for
1558 for( int i = 0; i < ims.planes(); ++i )
1559 {
1560 ims.image( i ) *= m_mask;
1561 }
1562 std::cerr << "done\n";
1563 }
1564
1565 if( m_preProcess_subradprof )
1566 {
1567 std::cerr << "subtracting radial profile . . .\n";
1568
1569#pragma omp parallel
1570 {
1571 eigenImageT rp;
1572
1573#pragma omp for
1574 for( int i = 0; i < ims.planes(); ++i )
1575 {
1576 Eigen::Map<Eigen::Array<realT, Eigen::Dynamic, Eigen::Dynamic>> imRef( ims.image( i ) );
1577 radprofim( rp, imRef, true );
1578 zeroNaNs( imRef );
1579 }
1580 }
1581
1582 std::cerr << "done\n";
1583
1584 if( m_maskFile != "" && m_preProcess_mask )
1585 {
1586 std::cerr << "Masking . . .\n";
1587#pragma omp parallel for
1588 for( int i = 0; i < ims.planes(); ++i )
1589 {
1590 ims.image( i ) *= m_mask;
1591 }
1592 std::cerr << "done\n";
1593 }
1594 }
1595
1596 if( m_preProcess_medianUSM_fwhm > 0 )
1597 {
1598 std::cerr << "applying median USM . . .\n";
1599
1600 eigenImageT medmask;
1601 medmask.resize( ims.rows(), ims.cols() );
1602 medmask.setConstant( 1 );
1603
1604 for( int z = 0; z < (int)( 0.5 * m_preProcess_medianUSM_fwhm ); ++z )
1605 {
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;
1610 }
1611
1612#pragma omp parallel
1613 {
1614 eigenImageT fim, im;
1615 fim.resize( ims.rows(), ims.cols() );
1616 im.resize( ims.rows(), ims.cols() );
1617
1618#pragma omp for
1619 for( int i = 0; i < ims.planes(); ++i )
1620 {
1621 im = ims.image( i );
1622 medianSmooth( fim, im, m_preProcess_medianUSM_fwhm );
1623 ims.image( i ) = ( im - fim ) * medmask;
1624 }
1625 }
1626
1627 if( m_maskFile != "" && m_preProcess_mask )
1628 {
1629 std::cerr << "Masking . . .\n";
1630#pragma omp parallel for
1631 for( int i = 0; i < ims.planes(); ++i )
1632 {
1633 ims.image( i ) *= m_mask;
1634 }
1635 }
1636 std::cerr << "done\n";
1637 }
1638
1639 if( m_preProcess_gaussUSM_fwhm > 0 )
1640 {
1641 std::cerr << "applying Gauss USM . . .\n";
1642 t_gaussusm_begin = sys::get_curr_time();
1643
1644#pragma omp parallel for
1645 for( int i = 0; i < ims.planes(); ++i )
1646 {
1647 eigenImageT fim, im;
1648 im = ims.image( i );
1649 filterImage( fim,
1650 im,
1651 gaussKernel<eigenImage<_realT>, 2>( m_preProcess_gaussUSM_fwhm ),
1652 0.5 * ( ims.cols() - 1 ) - m_preProcess_gaussUSM_fwhm * 4 );
1653 im = ( im - fim );
1654 ims.image( i ) = im;
1655 }
1656
1657 if( m_maskFile != "" && m_preProcess_mask )
1658 {
1659 std::cerr << "Masking . . .\n";
1660#pragma omp parallel for
1661 for( int i = 0; i < ims.planes(); ++i )
1662 {
1663 ims.image( i ) *= m_mask;
1664 }
1665 }
1666 t_gaussusm_end = sys::get_curr_time();
1667 std::cerr << "done\n";
1668 }
1669
1670 if( m_preProcess_azUSM_azW && m_preProcess_azUSM_radW )
1671 {
1672 ipc::ompLoopWatcher<> status( ims.planes(), std::cerr );
1673
1674 std::cerr << "applying azimuthal USM . . .\n";
1675 t_azusm_begin = sys::get_curr_time();
1676#pragma omp parallel for
1677 for( int i = 0; i < ims.planes(); ++i )
1678 {
1679 eigenImageT fim, im;
1680 im = ims.image( i );
1681 medianFilterImage( fim,
1682 im,
1683 azBoxKernel<eigenImage<realT>>( m_preProcess_azUSM_radW,
1684 m_preProcess_azUSM_azW,
1685 m_preProcess_azUSM_maxAz ) );
1686 im = ( im - fim );
1687 ims.image( i ) = im;
1688 status.incrementAndOutputStatus();
1689 }
1690
1691 if( m_maskFile != "" && m_preProcess_mask )
1692 {
1693 std::cerr << "Masking . . .\n";
1694#pragma omp parallel for
1695 for( int i = 0; i < ims.planes(); ++i )
1696 {
1697 ims.image( i ) *= m_mask;
1698 }
1699 }
1700
1701 t_azusm_end = sys::get_curr_time();
1702 std::cerr << "done (" << t_azusm_end - t_azusm_begin << " sec) \n";
1703 }
1704
1705 preProcess_meanSub( ims );
1706
1707 preProcess_pixelTSNorm( ims );
1708
1709 t_preproc_end = sys::get_curr_time();
1710
1711} // void HCIobservation<_realT>::preProcess()
1712
1713template <typename _realT>
1715{
1716 if( m_preProcess_meanSubMethod == HCI::meanSubMethod::none )
1717 {
1718 return;
1719 }
1720 else if( m_preProcess_meanSubMethod != HCI::meanSubMethod::meanImage &&
1721 m_preProcess_meanSubMethod != HCI::meanSubMethod::medianImage )
1722 {
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.";
1725 mxThrowException( err::invalidconfig, "HCIobservation::preProcess_meanSub", msg );
1726 }
1727
1728 eigenImageT mean;
1729
1730 if( m_preProcess_meanSubMethod == HCI::meanSubMethod::meanImage )
1731 {
1732 ims.mean( mean );
1733 }
1734 else if( m_preProcess_meanSubMethod == HCI::meanSubMethod::medianImage )
1735 {
1736 ims.median( mean );
1737 }
1738
1739#pragma omp parallel for
1740 for( int n = 0; n < ims.planes(); ++n )
1741 {
1742 ims.image( n ) -= mean;
1743
1744 if( m_maskFile != "" && m_preProcess_mask )
1745 {
1746 ims.image( n ) *= m_mask;
1747 }
1748 }
1749}
1750
1751template <typename _realT>
1753{
1754 if( m_preProcess_pixelTSNormMethod == HCI::pixelTSNormMethod::none )
1755 {
1756 return;
1757 }
1758
1759 if( m_preProcess_pixelTSNormMethod == HCI::pixelTSNormMethod::unknown )
1760 {
1761 mxThrowException( err::invalidconfig, "KlipReduction::preProcess_pixelTSNorm", "pixelTSNormMethod is unknown" );
1762 }
1763
1764 if( m_preProcess_pixelTSNormMethod == HCI::pixelTSNormMethod::rmsSigmaClipped )
1765 {
1766 mxThrowException( err::invalidconfig,
1767 "KlipReduction::preProcess_pixelTSNorm",
1768 "pixelTSNormMethod is rmsSigmaClipped, which is not implemented" );
1769 }
1770
1771 std::cerr << "normalizing pixels\n";
1772
1773#pragma omp parallel
1774 {
1775 std::vector<realT> pixs( ims.planes() );
1776
1777#pragma omp for
1778 for( int cc = 0; cc < ims.cols(); ++cc )
1779 {
1780 for( int rr = 0; rr < ims.rows(); ++rr )
1781 {
1782 if( m_maskFile != "" && m_preProcess_mask )
1783 {
1784 if( m_mask( rr, cc ) == 0 )
1785 {
1786 continue;
1787 }
1788 }
1789
1790 // We bother to load a vector in prep to add sigma clipping later.
1791 for( int pp = 0; pp < ims.planes(); ++pp )
1792 {
1793 pixs[pp] = ims.image( pp )( rr, cc );
1794 }
1795
1796 realT sd = sqrt( math::vectorVariance( pixs, 0 ) );
1797
1798 if( sd == 0 )
1799 continue;
1800
1801 for( int pp = 0; pp < ims.planes(); ++pp )
1802 {
1803 ims.image( pp )( rr, cc ) /= sd;
1804 }
1805 }
1806 }
1807 }
1808}
1809
1810template <typename _realT>
1812{
1813 std::ifstream fin;
1814 std::string str;
1815
1816 if( m_weightFile == "" )
1817 {
1818 mxError( "HCIobservation::readWeights", MXE_PARAMNOTSET, "m_weightFile not set" );
1819 return -1;
1820 }
1821
1822 // Read the weight file and load it into a map
1823 std::vector<std::string> wfileNames;
1824 std::vector<realT> imW;
1825
1826 if( ioutils::readColumns( m_weightFile, wfileNames, imW ) < 0 )
1827 return -1;
1828
1829 if( imW.size() < m_fileList.size() )
1830 {
1831 mxError( "HCIobservation::readWeights", MXE_SIZEERR, "not enough weights specified" );
1832 return -1;
1833 }
1834
1835 std::map<std::string, realT> weights;
1836 for( size_t i = 0; i < wfileNames.size(); ++i )
1837 weights[ioutils::pathFilename( wfileNames[i].c_str() )] = imW[i];
1838
1839 m_comboWeights.resize( m_fileList.size() );
1840
1841 realT wi;
1842 realT weightSum = 0;
1843 for( size_t i = 0; i < m_fileList.size(); ++i )
1844 {
1845 try
1846 {
1847 wi = weights.at( ioutils::pathFilename( m_fileList[i].c_str() ) );
1848 }
1849 catch( ... )
1850 {
1851 mxError( "HCIobservation::readWeights", MXE_NOTFOUND, "Weight for a file in m_fileList not found." );
1852 return -1;
1853 }
1854 m_comboWeights[i] = wi;
1855 weightSum += wi;
1856 }
1857
1858 // Finally normalize the weights
1859 for( size_t i = 0; i < m_comboWeights.size(); ++i )
1860 {
1861 m_comboWeights[i] /= weightSum;
1862 }
1863
1864 return 0;
1865} // int HCIobservation<_realT>::readWeights()
1866
1867template <typename _realT>
1869{
1870 if( m_combineMethod == HCI::noCombine )
1871 {
1872 return;
1873 }
1874
1875 t_combo_begin = sys::get_curr_time();
1876
1877 // Validate the combineMethod setting
1878 int method = HCI::medianCombine;
1879
1880 if( m_combineMethod == HCI::medianCombine )
1881 {
1882 method = HCI::medianCombine;
1883 }
1884 else if( m_combineMethod == HCI::meanCombine )
1885 {
1886 method = HCI::meanCombine;
1887 }
1888 else if( m_combineMethod == HCI::sigmaMeanCombine )
1889 {
1890 if( m_sigmaThreshold > 0 )
1891 {
1892 method = HCI::sigmaMeanCombine;
1893 }
1894 else
1895 {
1896 method = HCI::meanCombine;
1897 }
1898 }
1899 else if( m_combineMethod == HCI::debug )
1900 {
1901 method = HCI::debug;
1902 }
1903
1904 // Create and size temporary image for averaging
1905 eigenImageT tfinim;
1906
1907 m_finim.resize( m_psfsub[0].rows(), m_psfsub[0].cols(), m_psfsub.size() );
1908
1909 // Now cycle through each set of psf subtractions
1910 for( size_t n = 0; n < m_psfsub.size(); ++n )
1911 {
1912 if( method == HCI::medianCombine )
1913 {
1914 m_psfsub[n].median( tfinim );
1915 m_finim.image( n ) = tfinim;
1916 }
1917 else if( method == HCI::meanCombine )
1918 {
1919 if( m_comboWeights.size() == (size_t)m_Nims )
1920 {
1921 if( m_maskFile != "" )
1922 {
1923 m_psfsub[n].mean( tfinim, m_comboWeights, m_maskCube, m_minGoodFract );
1924 }
1925 else
1926 {
1927 m_psfsub[n].mean( tfinim, m_comboWeights );
1928 }
1929 }
1930 else
1931 {
1932 if( m_maskFile != "" )
1933 {
1934 m_psfsub[n].mean( tfinim, m_maskCube, m_minGoodFract );
1935 }
1936 else
1937 {
1938 m_psfsub[n].mean( tfinim );
1939 }
1940 }
1941 m_finim.image( n ) = tfinim;
1942 }
1943 else if( method == HCI::sigmaMeanCombine )
1944 {
1945 if( m_comboWeights.size() == (size_t)m_Nims )
1946 {
1947 if( m_maskFile != "" )
1948 {
1949 m_psfsub[n].sigmaMean( tfinim, m_comboWeights, m_maskCube, m_sigmaThreshold, m_minGoodFract );
1950 }
1951 else
1952 {
1953 m_psfsub[n].sigmaMean( tfinim, m_comboWeights, m_sigmaThreshold );
1954 }
1955 }
1956 else
1957 {
1958 if( m_maskFile != "" )
1959 {
1960 m_psfsub[n].sigmaMean( tfinim, m_maskCube, m_sigmaThreshold, m_minGoodFract );
1961 }
1962 else
1963 {
1964 m_psfsub[n].sigmaMean( tfinim, m_sigmaThreshold );
1965 }
1966 }
1967 m_finim.image( n ) = tfinim;
1968 }
1969 else if( method == HCI::debug )
1970 {
1971 m_finim.image( n ) = m_psfsub[n].image( 0 );
1972 }
1973 }
1974
1975 t_combo_end = sys::get_curr_time();
1976} // void HCIobservation<_realT>::combineFinim()
1977
1978template <typename _realT>
1980{
1981 if( m_preProcess_outputPrefix == "" )
1982 {
1983 return;
1984 }
1985
1986 std::string dir = ioutils::parentPath( m_preProcess_outputPrefix );
1988
1989 std::string fname;
1990
1992
1993 /** \todo Should add a HISTORY card here */
1994 char nstr[16];
1995 for( int i = 0; i < m_Nims; ++i )
1996 {
1997 snprintf( nstr, sizeof( nstr ), "%06d", i );
1998 fname = m_preProcess_outputPrefix + nstr + ".fits";
1999
2000 fits::fitsHeader fh = m_heads[i];
2001 stdFitsHeader( fh );
2002 ff.write( fname, m_tgtIms.image( i ).data(), m_Ncols, m_Nrows, 1, fh );
2003 }
2004} // void HCIobservation<_realT>::outputPreProcessed()
2005
2006template <typename _realT>
2008{
2009 head.append( "", fits::fitsCommentType(), "----------------------------------------" );
2010 head.append( "", fits::fitsCommentType(), "mx::HCIobservation parameters:" );
2011 head.append( "", fits::fitsCommentType(), "----------------------------------------" );
2012
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" );
2015
2016 head.append( "QFILE", ioutils::pathFilename( m_qualityFile.c_str() ), "quality file for thresholding" );
2017 head.append<realT>( "QTHRESH", m_qualityThreshold, "quality threshold" );
2018 head.append<int>( "NUMIMS", m_Nims, "number of images processed" );
2019
2020 head.append<int>( "IMSIZE", m_imSize, "image size after reading" );
2021
2022 head.append<std::string>( "COADMTHD", HCI::combineMethodStr( m_coaddCombineMethod ), "coadd combination method" );
2023 if( m_coaddCombineMethod != HCI::noCombine )
2024 {
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" );
2027 }
2028 else
2029 {
2030 head.append<int>( "COADIMNO", 0, "max number of images in each coadd" );
2031 head.append<realT>( "COADTIME", 0, "max time in each coadd" );
2032 }
2033
2034 head.append( "MASKFILE", m_maskFile, "mask file" );
2035
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" );
2039 head.append<realT>( "PREPROC AZUSM AZWIDTH",
2040 m_preProcess_azUSM_azW,
2041 "pre-process azimuthal USM azimuthal width [pixels]" );
2042 head.append<realT>( "PREPROC AZUSM MAXAZ",
2043 m_preProcess_azUSM_maxAz,
2044 "pre-process azimuthal USM maximum azimuthal width [degrees]" );
2045 head.append<realT>( "PREPROC AZUSM RADWIDTH",
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" );
2056}
2057
2058template <typename _realT>
2060{
2061 std::string fname = m_finimName;
2062
2063 if( m_outputDir != "" )
2064 {
2065 mkdir( m_outputDir.c_str(), S_IRWXU | S_IRWXG | S_IROTH | S_IXOTH );
2066
2067 fname = m_outputDir + "/" + fname;
2068 }
2069
2070 if( !m_exactFinimName )
2071 {
2072 fname = ioutils::getSequentialFilename( fname, ".fits" );
2073 }
2074
2075 fits::fitsHeader head;
2076
2077 // Add HCIobservation standard header:
2078 stdFitsHeader( head );
2079
2080 // Now add the final combination details:
2081 head.append<std::string>( "COMBINATION METHOD", HCI::combineMethodStr( m_combineMethod ), "combination method" );
2082
2083 if( m_weightFile != "" )
2084 head.append( "WEIGHT FILE", m_weightFile, "file containing weights for combination" );
2085
2086 if( m_combineMethod == HCI::sigmaMeanCombine )
2087 head.append<realT>( "SIGMA THRESHOLD", m_sigmaThreshold, "threshold for sigma clipping" );
2088
2089 head.append<realT>( "MIN FOOD FRACTION", m_minGoodFract, "minimum good fraction for inclusion" );
2090 if( addHead )
2091 {
2092 head.append( *addHead );
2093 }
2094
2095 fits::fitsHeaderGitStatus( head, "mxlib_uncomp", MXLIB_UNCOMP_CURRENT_SHA1, MXLIB_UNCOMP_REPO_MODIFIED );
2096
2098
2099 f.write( fname, m_finim, head );
2100
2101 std::cerr << "Final image written to: " << fname << "\n";
2102} // void HCIobservation<_realT>::writeFinim(fits::fitsHeader * addHead)
2103
2104template <typename _realT>
2106{
2107
2108 std::string fname;
2109
2110 fits::fitsHeader head;
2111
2112 // Add the HCIobservation standard fits header
2113 stdFitsHeader( head );
2114
2115 if( addHead )
2116 {
2117 head.append( *addHead );
2118 }
2119
2120 fits::fitsHeaderGitStatus( head, "mxlib_uncomp", MXLIB_UNCOMP_CURRENT_SHA1, MXLIB_UNCOMP_REPO_MODIFIED );
2121
2123
2124 std::ofstream wout;
2125
2126 if( m_comboWeights.size() > 0 )
2127 {
2128 fname = m_PSFSubPrefix + "weights.dat";
2129 if( m_outputDir != "" )
2130 {
2131 mkdir( m_outputDir.c_str(), S_IRWXU | S_IRWXG | S_IROTH | S_IXOTH );
2132 fname = m_outputDir + "/" + fname;
2133 }
2134 wout.open( fname );
2135 std::cerr << "writing comboWeights: " << fname << "\n";
2136 }
2137
2138 char num[256];
2139 for( size_t n = 0; n < m_psfsub.size(); ++n )
2140 {
2141 for( int p = 0; p < m_psfsub[n].planes(); ++p )
2142 {
2143 snprintf( num, 256, "_%03zu_%05d.fits", n, p );
2144 fname = m_PSFSubPrefix + num;
2145
2146 if( m_outputDir != "" )
2147 {
2148 mkdir( m_outputDir.c_str(), S_IRWXU | S_IRWXG | S_IROTH | S_IXOTH );
2149 fname = m_outputDir + "/" + fname;
2150 }
2151
2152 fits::fitsHeader h = head;
2153
2154 h.append( m_heads[p] );
2155
2156 f.write( fname, m_psfsub[n].image( p ).data(), m_psfsub[n].rows(), m_psfsub[n].cols(), 1, h );
2157
2158 if( m_comboWeights.size() > 0 && n == 0 )
2159 {
2160 wout << fname << " " << m_comboWeights[p] << "\n";
2161 }
2162 }
2163 }
2164
2165 if( m_comboWeights.size() > 0 )
2166 {
2167 wout.close();
2168 }
2169} // void HCIobservation<_realT>::outputPSFSub(fits::fitsHeader * addHead)
2170
2171template <typename _realT>
2172int HCIobservation<_realT>::readPSFSub( const std::string &dir,
2173 const std::string &prefix,
2174 const std::string &ext,
2175 size_t nReductions )
2176{
2177
2178 m_psfsub.resize( nReductions );
2179
2180 // Load first file to condigure based on its header.
2181 std::vector<std::string> flist = ioutils::getFileNames( dir, prefix, "000", ext );
2185
2186 ff.read( im, fh, flist[0] );
2187
2188 if( fh.count( "FDELFRNT" ) == 0 )
2189 {
2190 mxError( "KLIPReduction", MXE_PARAMNOTSET, "FDELFRNT not found in FITS header." );
2191 return -1;
2192 }
2193 m_deleteFront = fh["FDELFRNT"].value<int>();
2194 std::cerr << "deleteFront: " << m_deleteFront << "\n";
2195
2196 if( fh.count( "FDELBACK" ) == 0 )
2197 {
2198 mxError( "KLIPReduction", MXE_PARAMNOTSET, "FDELBACK not found in FITS header." );
2199 return -1;
2200 }
2201 m_deleteBack = fh["FDELBACK"].value<int>();
2202 std::cerr << "deleteBack: " << m_deleteBack << "\n";
2203
2204 if( fh.count( "QFILE" ) == 0 )
2205 {
2206 mxError( "KLIPReduction", MXE_PARAMNOTSET, "QFILE not found in FITS header." );
2207 return -1;
2208 }
2209 m_qualityFile = fh["QFILE"].String();
2210 std::cerr << "qualityFile: " << m_qualityFile << "\n";
2211
2212 if( fh.count( "QTHRESH" ) == 0 )
2213 {
2214 mxError( "KLIPReduction", MXE_PARAMNOTSET, "QTHRESH not found in FITS header." );
2215 return -1;
2216 }
2217 m_qualityThreshold = fh["QTHRESH"].value<realT>();
2218 std::cerr << "qualityThreshold: " << m_qualityThreshold << "\n";
2219
2220 if( fh.count( "COADMTHD" ) == 0 )
2221 {
2222 mxError( "KLIPReduction", MXE_PARAMNOTSET, "COADMTHD not found in FITS header." );
2223 return -1;
2224 }
2225 m_coaddCombineMethod = HCI::combineMethodFmStr( fh["COADMTHD"].String() );
2226 std::cerr << "coaddCombineMethod: " << m_coaddCombineMethod << "\n";
2227
2228 if( fh.count( "COADIMNO" ) != 0 )
2229 {
2230 m_coaddMaxImno = fh["COADIMNO"].value<int>();
2231 std::cerr << "coaddMaxImno: " << m_coaddMaxImno << "\n";
2232 }
2233
2234 if( fh.count( "COADTIME" ) != 0 )
2235 {
2236 m_coaddMaxImno = fh["COADTIME"].value<realT>();
2237 std::cerr << "coaddMaxtime: " << m_coaddMaxTime << "\n";
2238 }
2239
2240 if( m_maskFile == "" )
2241 {
2242 if( fh.count( "MASKFILE" ) == 0 )
2243 {
2244 mxError( "KLIPReduction",
2245 MXE_PARAMNOTSET,
2246 "MASKFILE not found in FITS header and not set in configuration." );
2247 return -1;
2248 }
2249 m_maskFile = fh["MASKFILE"].String();
2250 }
2251 std::cerr << "maskFile: " << m_maskFile << "\n";
2252
2253 if( fh.count( "PPBEFORE" ) == 0 )
2254 {
2255 mxError( "KLIPReduction", MXE_PARAMNOTSET, "PPBEFORE not found in FITS header." );
2256 return -1;
2257 }
2258 m_preProcess_beforeCoadd = fh["PPBEFORE"].value<int>();
2259 std::cerr << "preProcess_beforeCoadd: " << m_preProcess_beforeCoadd << "\n";
2260
2261 if( fh.count( "PPMASK" ) == 0 )
2262 {
2263 mxError( "KLIPReduction", MXE_PARAMNOTSET, "PPMASK not found in FITS header." );
2264 return -1;
2265 }
2266 m_preProcess_mask = fh["PPMASK"].value<int>();
2267 std::cerr << "preProcess_mask: " << m_preProcess_mask << "\n";
2268
2269 if( fh.count( "PPSUBRAD" ) == 0 )
2270 {
2271 mxError( "KLIPReduction", MXE_PARAMNOTSET, "PPSUBRAD not found in FITS header." );
2272 return -1;
2273 }
2274 m_preProcess_subradprof = fh["PPSUBRAD"].value<int>();
2275 std::cerr << "preProcess_subradprof: " << m_preProcess_subradprof << "\n";
2276
2277 if( fh.count( "PPAUSMAW" ) == 0 )
2278 {
2279 mxError( "KLIPReduction", MXE_PARAMNOTSET, "PPAUSMAW not found in FITS header." );
2280 return -1;
2281 }
2282 m_preProcess_azUSM_azW = fh["PPAUSMAW"].value<realT>();
2283 std::cerr << "preProcess_azUSM_azW: " << m_preProcess_azUSM_azW << "\n";
2284
2285 if( fh.count( "PPAUSMRW" ) == 0 )
2286 {
2287 mxError( "KLIPReduction", MXE_PARAMNOTSET, "PPAUSMRW not found in FITS header." );
2288 return -1;
2289 }
2290 m_preProcess_azUSM_radW = fh["PPAUSMRW"].value<realT>();
2291 std::cerr << "preProcess_azUSM_radW: " << m_preProcess_azUSM_radW << "\n";
2292
2293 if( fh.count( "PPGUSMFW" ) == 0 )
2294 {
2295 mxError( "KLIPReduction", MXE_PARAMNOTSET, "PPGUSMFW not found in FITS header." );
2296 return -1;
2297 }
2298 m_preProcess_gaussUSM_fwhm = fh["PPGUSMFW"].value<realT>();
2299 std::cerr << "preProcess_gaussUSM_fwhm: " << m_preProcess_gaussUSM_fwhm << "\n";
2300
2301 fits::fitsHeader head;
2302
2303 if( m_MJDKeyword != "" )
2304 head.append( m_MJDKeyword );
2305
2306 for( size_t i = 0; i < m_keywords.size(); ++i )
2307 {
2308 head.append( m_keywords[i] );
2309 }
2310
2311 for( size_t n = 0; n < nReductions; ++n )
2312 {
2313 char nstr[5];
2314 int nwr = snprintf( nstr, sizeof( nstr ), "%03zu", n );
2315 if( nwr < 0 || n >= sizeof( nstr ) )
2316 {
2317 std::cerr << "possibly bad formatting in filename\n";
2318 }
2319
2320 std::string nprefix = prefix + "_" + nstr + "_";
2321 load_fileList( dir, nprefix, ext );
2322
2323 if( m_fileList.size() == 0 )
2324 {
2325 mxError( "HCIobservation",
2326 MXE_FILENOTFOUND,
2327 "The m_fileList has 0 length, there are no files to be read." );
2328 return -1;
2329 }
2330
2331 Eigen::Array<realT, Eigen::Dynamic, Eigen::Dynamic> im;
2332
2333 fits::fitsFile<realT> f( m_fileList[0] );
2334
2335 fits::fitsHeader fh = head;
2336 f.read( im, fh );
2337
2338 // We set imSize to match the first image, but we make it a square.
2339 if( m_imSize == 0 )
2340 {
2341 m_imSize = im.rows();
2342 if( m_imSize > im.cols() )
2343 {
2344 m_imSize = im.cols();
2345 }
2346 }
2347 else
2348 {
2349 // Now make sure we don't read too much.
2350 if( m_imSize > im.rows() )
2351 {
2352 m_imSize = im.rows();
2353 }
2354 if( m_imSize > im.cols() )
2355 {
2356 m_imSize = im.cols();
2357 }
2358 }
2359
2360 // the +0.1 is just to make sure we don't have a problem with precision (we shouldn't)/
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 ),
2363 m_imSize,
2364 m_imSize );
2365 im.resize( m_imSize, m_imSize );
2366
2367 std::cerr << "image size is " << m_imSize << "\n";
2368
2369 if( n > 0 )
2370 {
2371 if( m_fileList.size() != (size_t)m_Nims )
2372 {
2373 mxError( "HCIobservation", MXE_INVALIDARG, "Different number of images in reductions." );
2374 return -1;
2375 }
2376 if( m_Nrows != im.rows() )
2377 {
2378 mxError( "HCIobservation", MXE_INVALIDARG, "Different number of rows in reductions." );
2379 return -1;
2380 }
2381 if( m_Ncols != im.cols() )
2382 {
2383 mxError( "HCIobservation", MXE_INVALIDARG, "Different number of cols in reductions." );
2384 return -1;
2385 }
2386 }
2387 else
2388 {
2389 std::cerr << "found " << nReductions << " sets of " << m_fileList.size() << " " << im.rows() << " x "
2390 << im.cols() << " files\n";
2391 }
2392 m_Nims = m_fileList.size();
2393 m_Nrows = im.rows();
2394 m_Ncols = im.cols();
2395 m_Npix = im.rows() * im.cols();
2396
2397 m_psfsub[n].resize( m_Nrows, m_Ncols, m_Nims );
2398
2399 m_heads.clear(); // This is necessary to make sure heads.resize() copies head on a 2nd call
2400 m_heads.resize( m_fileList.size(), head );
2401
2402 t_load_begin = sys::get_curr_time();
2403
2404 f.read( m_psfsub[n].data(), m_heads, m_fileList );
2405
2406 f.setReadSize();
2407
2408 if( m_MJDKeyword != "" )
2409 {
2410 m_imageMJD.resize( m_heads.size() );
2411
2412 if( m_MJDisISO8601 )
2413 {
2414 for( size_t i = 0; i < m_imageMJD.size(); ++i )
2415 {
2416 m_imageMJD[i] = sys::ISO8601date2mjd( m_heads[i][m_MJDKeyword].String() );
2417 }
2418 }
2419 else
2420 {
2421 for( size_t i = 0; i < m_imageMJD.size(); ++i )
2422 {
2423 m_imageMJD[i] = m_heads[i][m_MJDKeyword].template value<realT>() * m_MJDUnits;
2424 }
2425 }
2426 }
2427
2428 t_load_end = sys::get_curr_time();
2429
2430 for( size_t n = 0; n < m_psfsub.size(); ++n )
2431 {
2432 zeroNaNCube( m_psfsub[n] );
2433 }
2434 }
2435
2436 if( m_weightFile != "" )
2437 {
2438 std::vector<std::string> fn;
2439 ioutils::readColumns( m_weightFile, fn, m_comboWeights );
2440
2441 std::cerr << "read: " << m_weightFile << " (" << m_comboWeights.size() << ")\n";
2442 }
2443
2444 /*** Now do the post-read actions ***/
2445 if( postReadFiles() < 0 )
2446 {
2447 return -1;
2448 }
2449
2450 readMask();
2451
2452 m_filesRead = true;
2453
2454 return 0;
2455}
2456
2457///@}
2458
2459extern template class HCIobservation<float>;
2460extern template class HCIobservation<double>;
2461
2462} // namespace improc
2463} // namespace mx
2464
2465#endif //__HCIobservation_hpp__
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.
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.
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 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.
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.
Definition fileUtils.cpp:52
std::string parentPath(const std::string &fname)
Get the parent path from a filename.
Definition fileUtils.cpp:77
std::string pathFilename(const std::string &fname)
Get the base filename.
Definition fileUtils.cpp:71
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.
Definition timeUtils.hpp:90
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.
Image interpolation and transformation.
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 mxlib c++ namespace.
Definition mxError.hpp:106
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.
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.
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.