Citlali
Loading...
Searching...
No Matches
engine.h
Go to the documentation of this file.
1#pragma once
2
3#include "sys/types.h"
4#include "sys/sysinfo.h"
5
6#include <memory>
7#include <string>
8#include <vector>
9#include <omp.h>
10#include <fstream>
11#include <limits>
12
13#include <Eigen/Core>
14
15#include <citlali_config/config.h>
16#include <citlali_config/gitversion.h>
17#include <citlali_config/default_config.h>
18#include <kids/core/kidsdata.h>
19#include <kids/sweep/fitter.h>
20#include <kids/timestream/solver.h>
21#include <kids/toltec/toltec.h>
22#include <kidscpp_config/gitversion.h>
23#include <tula_config/gitversion.h>
24#include <tula/cli.h>
25#include <tula/config/core.h>
26#include <tula/config/flatconfig.h>
27#include <tula/config/yamlconfig.h>
28#include <tula/enum.h>
29#include <tula/filesystem.h>
30#include <tula/formatter/container.h>
31#include <tula/formatter/enum.h>
32#include <tula/grppi.h>
33#include <tula/logging.h>
34#include <tula/switch_invoke.h>
35
41
47
54
57
60
66
70
72 // interpolate over gaps in timestreams
74 // create reduction subdirectories
76
77 // run or skip tod processing
78 bool run_tod;
79
80 // output timestreams
82
83 // controls for mapmaking
88
89 // run source finding
91};
92
114
116 // source name
118
119 // beammap source position
121
122 // fluxes and errs
123 std::map<std::string, double> beammap_fluxes_mJy_beam, beammap_err_mJy_beam;
124 std::map<std::string, double> beammap_fluxes_MJy_Sr, beammap_err_MJy_Sr;
125
126 // maximum beammap iterations
128
129 // beammap tolerance
131
132 // subtract reference detector
134
135 // beammap reference detector
137
138 // derotate fitted detectors
140
141 // iteration to write out ptcdata
143
144 // upper and lower limits of psd for sensitivity calc
145 Eigen::VectorXd sens_psd_limits_Hz;
146
147 // limits on fwhm, sig2noise, and distance from center for flagging
148 std::map<std::string, double> lower_fwhm_arcsec, upper_fwhm_arcsec, lower_sig2noise,
150
151 // limits on sensitivity for flagging
153};
154
155class Engine: public reduControls, public reduClasses, public beammapControls {
156public:
157 // type for missing/invalid keys
158 using key_vec_t = std::vector<std::vector<std::string>>;
159
160 // get logger
161 std::shared_ptr<spdlog::logger> logger = spdlog::get("citlali_logger");
162
163 // for timing
164 Eigen::VectorXd t_common;
165 std::vector<Eigen::VectorXi> masks;
166 std::vector<Eigen::VectorXd> nw_times;
167
168 // date/time of each obs
169 std::vector<std::string> date_obs;
170
171 // add extra output for debugging
173
174 // time gaps
175 std::map<std::string,int> gaps;
176
177 // output directory and optional sub directory name
179
180 // reduction directory number
182
183 // obsnum and coadded directory names
185
186 // tod output file name
187 std::map<std::string, std::string> tod_filename;
188
189 // vectors to hold missing/invalid keys
191
192 // number of threads
194
195 // parallel execution policy
196 std::string parallel_policy;
197
198 // number of scans completed
200
201 // manual offsets for nws and hwp
202 std::map<std::string,double> interface_sync_offset;
203
204 // vectors for tod alignment offsets
205 std::vector<Eigen::Index> start_indices, end_indices;
206
207 // indices for hwpr alignment offsets
209
210 // xs, rs, is, qs
211 std::string tod_type;
212
213 // reduction type (science, pointing, beammap)
214 std::string redu_type;
215
216 // obsnum
217 std::string obsnum;
218
219 // write filtered maps as they complete
221
222 // rtc or ptc types
224
225 // map grouping and algorithm
227
228 // number of maps
230
231 // mapping from index in map vector to array index
233
234 // mapping from index in map vector to array index
235 Eigen::VectorXI maps_to_stokes;
236
237 // current fruit loops iteration
239
240 // manual pointing offsets
241 std::map<std::string, Eigen::VectorXd> pointing_offsets_arcsec;
242 // modified julian dates of pointing offsets
244
245 // map output files
246 std::vector<fitsIO<file_type_enum::write_fits, CCfits::ExtHDU*>> fits_io_vec, noise_fits_io_vec;
247 std::vector<fitsIO<file_type_enum::write_fits, CCfits::ExtHDU*>> filtered_fits_io_vec, filtered_noise_fits_io_vec;
248
249 // coadded map output files
250 std::vector<fitsIO<file_type_enum::write_fits, CCfits::ExtHDU*>> coadd_fits_io_vec, coadd_noise_fits_io_vec;
251 std::vector<fitsIO<file_type_enum::write_fits, CCfits::ExtHDU*>> filtered_coadd_fits_io_vec, filtered_coadd_noise_fits_io_vec;
252
253 // per obsnum setup common to all redu types
254 void obsnum_setup();
255
256 // get RTC config options
257 template<typename CT>
258 void get_rtc_config(CT &);
259
260 // get PTC config options
261 template<typename CT>
262 void get_ptc_config(CT &);
263
264 // get timestream config options
265 template<typename CT>
266 void get_timestream_config(CT &);
267
268 // get beammap config options
269 template<typename CT>
270 void get_beammap_config(CT &);
271
272 // get mapmaking config options
273 template<typename CT>
274 void get_mapmaking_config(CT &);
275
276 // get map filtering config options
277 template<typename CT>
278 void get_map_filter_config(CT &);
279
280 // get all non-input config options and call other config functions
281 template<typename CT>
282 void get_citlali_config(CT &);
283
284 // get source fluxes (beammap only)
285 template<typename CT>
286 void get_photometry_config(CT &);
287
288 // get pointing offsets
289 template<typename CT>
290 void get_astrometry_config(CT &);
291
292 // create fits files (does not populate them)
294
295 // add FITS header values to tod files
296 void add_tod_header();
297
298 // create tod files (does not populate them)
299 template <engine_utils::toltecIO::ProdType prod_t>
300 void create_tod_files();
301
302 // output obs summary at command line
303 void cli_summary();
304
305 // write time chunk summary (verbose mode)
306 template <TCDataKind tc_t>
308
309 // write map summary (verbose mode)
310 template <typename map_buffer_t>
311 void write_map_summary(map_buffer_t &);
312
313 // create filenames
316 auto setup_filenames(std::string dir_name);
317
318 // create variable names for maps, psds, and hists
319 auto get_map_name(int);
320
321 // add primary header to FITS files
322 template <typename fits_io_type, class map_buffer_t>
323 void add_phdu(fits_io_type &, map_buffer_t &, Eigen::Index);
324
325 // add maps to FITS files and output them
326 template <typename fits_io_type, class map_buffer_t>
327 void write_maps(fits_io_type &, fits_io_type &, map_buffer_t &, Eigen::Index);
328
329 // write map psds
330 template <mapmaking::MapType map_t, class map_buffer_t>
331 void write_psd(map_buffer_t &, std::string);
332
333 // write map histograms
334 template <mapmaking::MapType map_t, class map_buffer_t>
335 void write_hist(map_buffer_t &, std::string);
336
337 // write stats netCDF4 file
338 void write_stats();
339
340 // run the wiener filter
341 template <mapmaking::MapType map_t, class map_buffer_t>
342 void run_wiener_filter(map_buffer_t &);
343
344 // find sources in the maps
345 template <mapmaking::MapType map_t, class map_buffer_t>
346 void find_sources(map_buffer_t &);
347
348 // write the sources to ecsv table
349 template <mapmaking::MapType map_t, class map_buffer_t>
350 void write_sources(map_buffer_t &, std::string);
351};
352
355 // get atm model
357
358 logger->info("using {} model for extinction correction",rtcproc.calibration.extinction_model);
359
360 // check tau (may be unnecessary now)
361 if (!telescope.sim_obs) {
362 Eigen::VectorXd tau_el(1);
363 // get mean elevation
364 tau_el << telescope.tel_data["TelElAct"].mean();
365 // get tau at mean elevation for each band
366 auto tau_freq = rtcproc.calibration.calc_tau(tau_el, telescope.tau_225_GHz);
367 // loop through and make sure average tau is not negative (implies wrong model)
368 for (auto const& [key, val] : tau_freq) {
369 if (val[0] < 0) {
370 logger->error("calculated mean {} tau {} < 0",toltec_io.array_name_map[key], val[0]);
371 std::exit(EXIT_FAILURE);
372 }
373 }
374 }
375 }
376 else {
378 }
379
380 // make sure there are matched fg's in apt if reducing in polarized mode
382 if ((calib.apt["fg"].array()==-1).all()) {
383 logger->error("no matched freq groups. cannot run in polarized mode");
384 std::exit(EXIT_FAILURE);
385 }
386 }
387
388 // setup kernel
389 if (rtcproc.run_kernel) {
391 }
392
393 // set despiker sample rate
395
396 // if filter is requested, make it here
399 }
400
401 // set map wcs crvals to source ra/dec
402 if (telescope.pixel_axes == "radec") {
403 omb.wcs.crval[0] = telescope.tel_header["Header.Source.Ra"](0)*RAD_TO_DEG;
404 omb.wcs.crval[1] = telescope.tel_header["Header.Source.Dec"](0)*RAD_TO_DEG;
405
406 if (run_coadd) {
407 cmb.wcs.crval[0] = telescope.tel_header["Header.Source.Ra"](0)*RAD_TO_DEG;
408 cmb.wcs.crval[1] = telescope.tel_header["Header.Source.Dec"](0)*RAD_TO_DEG;
409 }
410 }
411
412 // set map wcs crvals to source l/b
413 else if (telescope.pixel_axes == "galactic") {
414 omb.wcs.crval[0] = telescope.tel_header["Header.Source.L"](0)*RAD_TO_DEG;
415 omb.wcs.crval[1] = telescope.tel_header["Header.Source.B"](0)*RAD_TO_DEG;
416
417 if (run_coadd) {
418 cmb.wcs.crval[0] = telescope.tel_header["Header.Source.L"](0)*RAD_TO_DEG;
419 cmb.wcs.crval[1] = telescope.tel_header["Header.Source.B"](0)*RAD_TO_DEG;
420 }
421 }
422
423 // create timestream files
424 if (run_tod_output) {
425 // create tod output subdirectory if requested
426 if (tod_output_subdir_name!="null") {
427 fs::create_directories(obsnum_dir_name + "raw/" + tod_output_subdir_name);
428 }
429 // make rtc tod output file
430 if (tod_output_type == "rtc" || tod_output_type == "both") {
431 create_tod_files<engine_utils::toltecIO::rtc_timestream>();
432 }
433 // make ptc tod output file
434 if (tod_output_type == "ptc" || tod_output_type == "both") {
435 create_tod_files<engine_utils::toltecIO::ptc_timestream>();
436 }
437 }
438 // don't calculate any eigenvalues
439 else if (!diagnostics.write_evals) {
441 }
442
443 // tod output mode require sequential policy so set explicitly
444 if (run_tod_output) {
445 logger->warn("tod output mode require sequential policy. "
446 "parallelization will be disabled for some stages.");
447 parallel_policy = "seq";
448 }
449
450 // output basic info for obs reduction to command line
451 cli_summary();
452
453 // set up per-det stats file values
454 for (const auto &stat: diagnostics.det_stats_header) {
455 diagnostics.stats[stat].setZero(calib.n_dets, telescope.scan_indices.cols());
456 }
457 // set up per-group stats file values
458 for (const auto &stat: diagnostics.grp_stats_header) {
459 diagnostics.stats[stat].setZero(calib.n_arrays, telescope.scan_indices.cols());
460 }
461 // clear stored eigenvalues
462 std::map<Eigen::Index, std::vector<std::vector<Eigen::VectorXd>>>().swap(diagnostics.evals);
463}
464
465template<typename CT>
466void Engine::get_rtc_config(CT &config) {
467 logger->info("getting rtc config options");
468 // get rtcproc config
470
471 // offset inner chunks
474 }
475 // otherwise start at zero
476 else {
478 }
479
480 // ignore hwpr?
482 std::tuple{"timestream","polarimetry", "ignore_hwpr"});
483}
484
485template<typename CT>
486void Engine::get_ptc_config(CT &config) {
487 logger->info("getting ptc config options");
488 // get ptcproc config
490
491 // copy tod output bool for eigenvalues
494}
495
496template<typename CT>
498 logger->info("getting timestream config options");
499 // run tod processing
501 std::tuple{"timestream","enabled"});
502 // tod type (xs, rs, is, qs)
504 std::tuple{"timestream","type"});
505
506 // run rtc or ptc tod output?
507 bool run_tod_output_rtc, run_tod_output_ptc;
508 // output rtc
509 get_config_value(config, run_tod_output_rtc, missing_keys, invalid_keys,
510 std::tuple{"timestream","raw_time_chunk","output","enabled"});
511 // output ptc
512 get_config_value(config, run_tod_output_ptc, missing_keys, invalid_keys,
513 std::tuple{"timestream","processed_time_chunk","output","enabled"});
514 // set tod output to false by default
515 run_tod_output = false;
516
517 // check if rtc output is requested
518 if (run_tod_output_rtc) {
519 run_tod_output = true;
520 tod_output_type = "rtc";
521 }
522 // if ptc output is requested
523 if (run_tod_output_ptc) {
524 // check if rtc output was requested
525 if (run_tod_output == true) {
526 tod_output_type = "both";
527 }
528 // else just output ptc
529 else {
530 run_tod_output = true;
531 tod_output_type = "ptc";
532 }
533 }
534
535 // tod subdirectory name
537 std::tuple{"timestream","output", "subdir_name"});
538 // write eigenvalues to stats file
540 std::tuple{"timestream","output", "stats","eigenvalues"});
541 // get time chunk size
543 std::tuple{"timestream","chunking", "chunk_mode"});
544 // get time chunk size
546 std::tuple{"timestream","chunking", "value"});
547 // force chunking?
549 std::tuple{"timestream","chunking", "force_chunking"});
550
551 /* get raw time chunk config */
552 get_rtc_config(config);
553
554 /* get processed time chunk config */
555 get_ptc_config(config);
556}
557
558template<typename CT>
560 logger->info("getting mapmaking config options");
561 // enable mapmaking?
563 std::tuple{"mapmaking","enabled"});
564 // map grouping
566 std::tuple{"mapmaking","grouping"},{"auto","array","nw","detector","fg"});
567
568 // polarization is disabled for detector grouping
569 if (rtcproc.run_polarization && ((redu_type=="beammap" && map_grouping=="auto") || map_grouping=="detector")) {
570 logger->error("Detector grouping reductions do not currently support polarimetry mode");
571 std::exit(EXIT_FAILURE);
572 }
573
574 // set rtcproc map_grouping
576
577 // map_method
579 std::tuple{"mapmaking","method"},{"naive","jinc","maximum_likelihood"});
580
581 // map reference frame (radec, altaz, galactic)
583 std::tuple{"mapmaking","pixel_axes"},{"radec","altaz", "galactic"});
584
585 // get config for omb
586 logger->info("getting omb config options");
588
589 // run coaddition?
591 std::tuple{"coadd","enabled"});
592
593 // re-run to get config for cmb
594 if (run_coadd) {
595 logger->info("getting cmb config options");
597 }
598
599 // if flux calibration is not enabled, use tod type units (xs, rs, is, or qs)
600 if (!rtcproc.run_calibrate) {
603 }
604
605 // set parallelization for psd filter ffts (maintained with tod output/verbose mode)
609
610 if (map_method=="jinc") {
611 // maximum radius for jinc filter
613 std::tuple{"mapmaking","jinc_filter","r_max"});
614 // get jinc filter shape params
615 for (auto const& [arr_index, arr_name] : toltec_io.array_name_map) {
616 auto jinc_shape_vec = config.template get_typed<std::vector<double>>(std::tuple{"mapmaking","jinc_filter","shape_params",arr_name});
617 jinc_mm.shape_params[arr_index] = Eigen::Map<Eigen::VectorXd>(jinc_shape_vec.data(),jinc_shape_vec.size());
618 }
619
620 if (jinc_mm.mode=="matrix") {
621 // allocate jinc matrix
623 }
624 else if (jinc_mm.mode=="splines") {
625 // precompute jinc spline
627 }
628 }
629
630 else if (map_method=="maximum_likelihood") {
632 std::tuple{"mapmaking","maximum_likelihood","tolerance"});
634 std::tuple{"mapmaking","maximum_likelihood","max_iterations"});
635 }
636
637 // make noise maps?
639 std::tuple{"noise_maps","enabled"});
640 if (run_noise) {
641 // number of noise maps
643 std::tuple{"noise_maps","n_noise_maps"},{},{0},{});
644 // randomize noise maps on detector as well as time chunk
646 std::tuple{"noise_maps","randomize_dets"});
647
648 if (run_coadd) {
649 // copy omb number of noise maps to cmb
651 // copy randomize_dets to cmb
653 }
654 }
655 // otherwise set number of noise maps to zero
656 else {
657 omb.n_noise = 0;
658 cmb.n_noise = 0;
659 }
660
661 // set mapmaker polarization
664}
665
666template<typename CT>
668 logger->info("getting beammap config options");
669 // max beammap iteration
671 std::tuple{"beammap","iter_max"});
672 // beammap iteration tolerance
674 std::tuple{"beammap","iter_tolerance"});
675 // beammap reference detector
677 std::tuple{"beammap","reference_det"});
678 // subtract reference detector?
680 std::tuple{"beammap","subtract_reference_det"});
681 // derotate apt?
683 std::tuple{"beammap","derotate"});
684
685 // lower fwhm limit
686 auto lower_fwhm_arcsec_vec = config.template get_typed<std::vector<double>>(std::tuple{"beammap","flagging","array_lower_fwhm_arcsec"});
687 // upper fwhm limit
688 auto upper_fwhm_arcsec_vec = config.template get_typed<std::vector<double>>(std::tuple{"beammap","flagging","array_upper_fwhm_arcsec"});
689 // lower signal-to-noise limit
690 auto lower_sig2noise_vec = config.template get_typed<std::vector<double>>(std::tuple{"beammap","flagging","array_lower_sig2noise"});
691 // upper signal-to-noise limit
692 auto upper_sig2noise_vec = config.template get_typed<std::vector<double>>(std::tuple{"beammap","flagging","array_upper_sig2noise"});
693 // maximum allowed distance limit
694 auto max_dist_arcsec_vec = config.template get_typed<std::vector<double>>(std::tuple{"beammap","flagging","array_max_dist_arcsec"});
695
696 // add params to respective array values
697 Eigen::Index i = 0;
698 for (auto const& [arr_index, arr_name] : toltec_io.array_name_map) {
699 // lower fwhm limit
700 lower_fwhm_arcsec[arr_name] = lower_fwhm_arcsec_vec[i];
701 // upper fwhm limit
702 upper_fwhm_arcsec[arr_name] = upper_fwhm_arcsec_vec[i];
703 // lower signal-to-noise limit
704 lower_sig2noise[arr_name] = lower_sig2noise_vec[i];
705 // upper signal-to-noise limit
706 upper_sig2noise[arr_name] = upper_sig2noise_vec[i];
707 // maximum allowed distance limit
708 max_dist_arcsec[arr_name] = max_dist_arcsec_vec[i];
709 i++;
710 }
711
712 // sensitivity factors
713 auto sens_factors_vec = config.template get_typed<std::vector<double>>(std::tuple{"beammap","flagging","sens_factors"});
714 lower_sens_factor = sens_factors_vec[0];
715 upper_sens_factor = sens_factors_vec[1];
716
717 // upper and lower frequencies over which to calculate sensitivity
718 sens_psd_limits_Hz.resize(2);
719 // get psd limits for sens from config
720 auto sens_psd_limits_Hz_vec = config.template get_typed<std::vector<double>>(std::tuple{"beammap","sens_psd_limits_Hz"});
721 // map sens limits back to Eigen vector
722 sens_psd_limits_Hz = (Eigen::Map<Eigen::VectorXd>(sens_psd_limits_Hz_vec.data(), sens_psd_limits_Hz_vec.size()));
723
724 // if no tolerance is specified, write out max iteration tod
725 if (run_tod_output) {
726 if (beammap_iter_tolerance <=0) {
728 }
729 // otherwise write out first iteration tod
730 else {
732 }
733 }
734}
735
736template<typename CT>
738 logger->info("getting map filtering config options");
739 // get wiener filter config options
741
742 // if in science mode, write filtered maps as they complete
743 if (redu_type=="science") {
745 }
746 // otherwise write at end
747 else {
749 }
750 // check if kernel is enabled
751 if (wiener_filter.template_type=="kernel") {
752 if (!rtcproc.run_kernel) {
753 logger->error("wiener filter kernel template requires kernel");
754 std::exit(EXIT_FAILURE);
755 }
756 // copy the map fitter
757 else {
759 }
760 }
761 // make sure noise maps were enabled
762 if (!run_noise && (!wiener_filter.run_lowpass && wiener_filter.filter_type=="wiener_filter")) {
763 logger->error("wiener filter requires noise maps");
764 std::exit(EXIT_FAILURE);
765 }
766
767 // set parallelization for ffts (maintained with tod output/verbose mode)
769}
770
771template<typename CT>
773 // get interface offsets
774 if (config.has(std::tuple{"interface_sync_offset"})) {
775 auto interface_node = config.get_node(std::tuple{"interface_sync_offset"});
776 // interface key names
777 std::vector<std::string> interface_keys = {
778 "toltec0",
779 "toltec1",
780 "toltec2",
781 "toltec3",
782 "toltec4",
783 "toltec5",
784 "toltec6",
785 "toltec7",
786 "toltec8",
787 "toltec9",
788 "toltec10",
789 "toltec11",
790 "toltec12",
791 "hwpr"
792 };
793 // loop through interfaces
794 for (Eigen::Index i=0; i<interface_node.size(); ++i) {
795 auto offset = config.template get_typed<double>(std::tuple{"interface_sync_offset",i, interface_keys[i]});
796 interface_sync_offset[interface_keys[i]] = offset;
797 }
798 }
799
800 // verbose mode?
802 std::tuple{"runtime","verbose"});
803 // output directory
805 std::tuple{"runtime","output_dir"});
806 // number of threads to use
808 std::tuple{"runtime","n_threads"});
809 // overall parallel policy
811 std::tuple{"runtime","parallel_policy"},{"seq","omp"});
812 // reduction type (science, pointing, beammap)
814 std::tuple{"runtime","reduction_type"},{"science","pointing","beammap"});
815 // create redu00, redu01... subdirectories
817 std::tuple{"runtime","use_subdir"});
818 // interp over gaps in align_timestream
820 std::tuple{"runtime","interp_over_gaps"});
821
822 /* get timestream config */
823 get_timestream_config(config);
824
825 /* get mapmaking config */
826 get_mapmaking_config(config);
827
828 // run map filter?
830 std::tuple{"post_processing","map_filtering","enabled"});
831
832 // run source finder?
834 std::tuple{"post_processing","source_finding","enabled"});
835
836 // map fitter options if in pointing or beammap mode or if map filtering or source finding are enabled
837 if (redu_type=="pointing" || redu_type=="beammap" || run_map_filter || run_source_finder) {
838 // size of region around found source to fit
840 std::tuple{"post_processing","source_fitting","bounding_box_arcsec"},{},{0});
841 // radius around center of map to find source within
843 std::tuple{"post_processing","source_fitting","fitting_radius_arcsec"});
844 // fit 2d gaussian rotation angle
846 std::tuple{"post_processing","source_fitting", "gauss_model","fit_rotation_angle"});
847
848 // convert bounding box and fitting region to pixels
851
852 // fitter flux and fwhm limits
853 map_fitter.flux_limits.resize(2);
854 map_fitter.fwhm_limits.resize(2);
855 for (Eigen::Index i=0; i<map_fitter.flux_limits.size(); ++i) {
856 // flux limit
857 map_fitter.flux_limits(i) = config.template get_typed<double>(std::tuple{"post_processing","source_fitting",
858 "gauss_model","amp_limit_factors",i});
859 // fwhm limit
860 map_fitter.fwhm_limits(i) = config.template get_typed<double>(std::tuple{"post_processing","source_fitting",
861 "gauss_model","fwhm_limit_factors",i});
862 }
863
864 // flux lower factor
865 if (map_fitter.flux_limits(0) > 0) {
867 }
868 // flux lower factor
869 if (map_fitter.flux_limits(1) > 0) {
871 }
872 // fwhm lower factor
873 if (map_fitter.fwhm_limits(0) > 0) {
875 }
876 // fwhm upper factor
877 if (map_fitter.fwhm_limits(1) > 0) {
879 }
880 }
881
882 /* get wiener filter config */
883 if (run_map_filter) {
884 // needs map fitter config
885 get_map_filter_config(config);
886 }
887
888 // get source finder config options
889 if (run_source_finder) {
890 // minimum found source sigma
892 std::tuple{"post_processing","source_finding","source_sigma"});
893 // window around source to exclude other sources
895 std::tuple{"post_processing","source_finding","source_window_arcsec"});
896 // search map, negative of map, or both
898 std::tuple{"post_processing","source_finding","mode"});
899
900 // convert source window to radians
902
903 if (run_coadd) {
904 // copy omb source sigma to cmb
906 // copy omb source_window_rad to cmb
908 // copy omb source_finder_mode to cmb
910 }
911 }
912
913 /* get beammap config */
914 if (redu_type=="beammap") {
915 // needs redu_type config
916 get_beammap_config(config);
917 }
918
919 // disable map related keys if map-making is disabled
920 if (!run_mapmaking) {
921 run_coadd = false;
922 run_noise = false;
923 run_map_filter = false;
924 run_source_finder = false;
925 // we don't need to do iterations if no maps are made
927 }
928}
929
930template<typename CT>
932 // beammap source name
934 std::tuple{"beammap_source","name"});
935 // beammap source ra
937 std::tuple{"beammap_source","ra_deg"});
938 // convert ra to radians
940
941 // beammap source dec
943 std::tuple{"beammap_source","dec_deg"});
944 // convert dec to radians
946
947 // number of fluxes
948 Eigen::Index n_fluxes = config.get_node(std::tuple{"beammap_source","fluxes"}).size();
949
950 // get source fluxes
951 for (Eigen::Index i=0; i<n_fluxes; ++i) {
952 auto array = config.get_str(std::tuple{"beammap_source","fluxes",i,"array_name"});
953 // source flux in mJy/beam
954 auto flux = config.template get_typed<double>(std::tuple{"beammap_source","fluxes",i,"value_mJy"});
955 // source flux uncertainty in mJy/beam
956 auto uncertainty_mJy = config.template get_typed<double>(std::tuple{"beammap_source","fluxes",i,"uncertainty_mJy"});
957
958 // copy flux and uncertainty
959 beammap_fluxes_mJy_beam[array] = flux;
960 beammap_err_mJy_beam[array] = uncertainty_mJy;
961 }
962}
963
964template<typename CT>
966 // check if config file has pointing_offsets
967 if (config.has("pointing_offsets")) {
968 std::vector<double> offset;
969 // get az offset
970 offset = config.template get_typed<std::vector<double>>(std::tuple{"pointing_offsets",0,"value_arcsec"});
971 pointing_offsets_arcsec["az"] = Eigen::Map<Eigen::VectorXd>(offset.data(),offset.size());
972 // get alt offset
973 offset = config.template get_typed<std::vector<double>>(std::tuple{"pointing_offsets",1,"value_arcsec"});
974 pointing_offsets_arcsec["alt"] = Eigen::Map<Eigen::VectorXd>(offset.data(),offset.size());
975
976 // get julian date of pointing offsets
977 try {
978 offset = config.template get_typed<std::vector<double>>(std::tuple{"pointing_offsets",2,"modified_julian_date"});
979 pointing_offsets_modified_julian_date = Eigen::Map<Eigen::VectorXd>(offset.data(),offset.size());
980 }
981 catch (...) {
983 }
984 }
985 else {
986 logger->error("pointing_offsets not found in config");
987 std::exit(EXIT_FAILURE);
988 }
989}
990
992 // clear fits vectors for each observation
993 fits_io_vec.clear();
994 noise_fits_io_vec.clear();
995 filtered_fits_io_vec.clear();
997
998 // loop through arrays
999 for (Eigen::Index i=0; i<calib.n_arrays; ++i) {
1000 // array index
1001 auto array = calib.arrays[i];
1002 // array name
1003 std::string array_name = toltec_io.array_name_map[array];
1004 // map filename
1008 // create fits_io class for current array file
1010 // append to fits_io vector
1011 fits_io_vec.push_back(std::move(fits_io));
1012
1013 // if noise maps are requested but coadding is not, populate noise fits vector
1014 if (!run_coadd) {
1015 if (run_noise) {
1016 // noise map filename
1020 // create fits_io class for current array file
1022 // append to fits_io vector
1023 noise_fits_io_vec.push_back(std::move(fits_io));
1024 }
1025
1026 // map filtering
1027 if (run_map_filter) {
1028 // filtered map filename
1031 redu_type, array_name,
1033 // create fits_io class for current array file
1035 // append to fits_io vector
1036 filtered_fits_io_vec.push_back(std::move(fits_io));
1037
1038 // filtered noise maps
1039 if (run_noise) {
1040 // filtered noise map filename
1043 array_name, obsnum, telescope.sim_obs);
1044 // create fits_io class for current array file
1046 // append to fits_io vector
1047 filtered_noise_fits_io_vec.push_back(std::move(fits_io));
1048 }
1049 }
1050 }
1051 }
1052}
1053
1055 // loop through viles
1056 for (const auto & [fkey, fval]: tod_filename) {
1057 netCDF::NcFile fo(fval, netCDF::NcFile::write);
1058
1059 // add unit conversions
1060 if (rtcproc.run_calibrate) {
1061 for (const auto &val: calib.arrays) {
1062 auto name = toltec_io.array_name_map[val];
1063 // conversion to uK
1064 auto fwhm = (std::get<0>(calib.array_fwhms[val]) + std::get<1>(calib.array_fwhms[val]))/2;
1065 auto mJy_beam_to_uK = engine_utils::mJy_beam_to_uK(1, toltec_io.array_freq_map[val], fwhm*ASEC_TO_RAD);
1066
1067 // beam area in steradians
1068 auto beam_area_rad = 2.*pi*pow(fwhm*FWHM_TO_STD*ASEC_TO_RAD,2);
1069 // get Jy/pixel
1070 auto mJy_beam_to_Jy_px = 1e-3/beam_area_rad*pow(omb.pixel_size_rad,2);
1071
1072 if (omb.sig_unit == "mJy/beam") {
1073 // conversion to mJy/beam
1074 add_netcdf_var(fo, "to_mJy_beam_"+name, 1);
1075 // conversion to MJy/sr
1076 add_netcdf_var(fo, "to_MJy_sr_"+name, 1/(calib.array_beam_areas[val]*MJY_SR_TO_mJY_ASEC));
1077 // conversion to uK
1078 add_netcdf_var(fo, "to_uK_"+name, mJy_beam_to_uK);
1079 // conversion to Jy/pixel
1080 add_netcdf_var(fo, "to_Jy_pixel_"+name, mJy_beam_to_Jy_px);
1081 }
1082 else if (omb.sig_unit == "MJy/sr") {
1083 // conversion to mJy/beam
1084 add_netcdf_var(fo, "to_mJy_beam_"+name, calib.array_beam_areas[val]*MJY_SR_TO_mJY_ASEC);
1085 // conversion to MJy/Sr
1086 add_netcdf_var(fo, "to_MJy_sr_"+name, 1);
1087 // conversion to uK
1088 add_netcdf_var(fo, "to_uK_"+name, calib.array_beam_areas[val]*MJY_SR_TO_mJY_ASEC*mJy_beam_to_uK);
1089 // conversion to Jy/pixel
1090 add_netcdf_var(fo, "to_Jy_pixel_"+name, calib.array_beam_areas[val]*MJY_SR_TO_mJY_ASEC*mJy_beam_to_Jy_px);
1091 }
1092 else if (omb.sig_unit == "uK") {
1093 // conversion to mJy/beam
1094 add_netcdf_var(fo, "to_mJy_beam_"+name, 1/mJy_beam_to_uK);
1095 // conversion to MJy/sr
1096 add_netcdf_var(fo, "to_MJy_sr_"+name, 1/mJy_beam_to_uK/(calib.array_beam_areas[val]*MJY_SR_TO_mJY_ASEC));
1097 // conversion to uK
1098 add_netcdf_var(fo, "to_uK_"+name, 1);
1099 // conversion to Jy/pixel
1100 add_netcdf_var(fo, "to_Jy_pixel_"+name, (1/mJy_beam_to_uK)*mJy_beam_to_Jy_px);
1101 }
1102 else if (omb.sig_unit == "Jy/pixel") {
1103 // conversion to mJy/beam
1104 add_netcdf_var(fo, "to_mJy_beam_"+name, 1/mJy_beam_to_Jy_px);
1105 // conversion to MJy/sr
1106 add_netcdf_var(fo, "to_MJy_sr_"+name, (1/mJy_beam_to_Jy_px)/(calib.array_beam_areas[val]*MJY_SR_TO_mJY_ASEC));
1107 // conversion to uK
1108 add_netcdf_var(fo, "to_uK_"+name, mJy_beam_to_uK/mJy_beam_to_Jy_px);
1109 // conversion to Jy/pixel
1110 add_netcdf_var(fo, "to_Jy_pixel_"+name, 1);
1111 }
1112 }
1113 }
1114
1115 // add date and time of obs
1116 add_netcdf_var<std::string>(fo, "DATEOBS0", date_obs.back());
1117
1118 // add source
1119 add_netcdf_var<std::string>(fo,"SOURCE",telescope.source_name);
1120
1121 // add source flux for beammaps
1122 if (redu_type == "beammap") {
1123 for (const auto &val: calib.arrays) {
1124 auto name = toltec_io.array_name_map[val];
1125 add_netcdf_var(fo, "HEADER.SOURCE.FLUX_MJYPERBEAM_"+name, beammap_fluxes_mJy_beam[name]);
1126 add_netcdf_var(fo, "HEADER.SOURCE.FLUX_MJYPERSR_"+name, beammap_fluxes_MJy_Sr[name]);
1127 }
1128 add_netcdf_var(fo, "BEAMMAP.ITER_TOLERANCE", beammap_iter_tolerance);
1129 add_netcdf_var(fo, "BEAMMAP.ITER_MAX", beammap_iter_max);
1130 add_netcdf_var(fo, "BEAMMAP.IS_DEROTATED", beammap_derotate);
1131
1132 // add reference detector information
1134 add_netcdf_var(fo, "BEAMMAP.REF_DET_INDEX", beammap_reference_det);
1135 add_netcdf_var(fo, "BEAMMAP.REF_X_T", calib.apt["x_t"](beammap_reference_det));
1136 add_netcdf_var(fo, "BEAMMAP.REF_Y_T", calib.apt["y_t"](beammap_reference_det));
1137 }
1138 else {
1139 add_netcdf_var(fo, "BEAMMAP.REF_DET_INDEX", -99);
1140 add_netcdf_var(fo, "BEAMMAP.REF_X_T", -99);
1141 add_netcdf_var(fo, "BEAMMAP.REF_Y_T", -99);
1142 }
1143 }
1144
1145 add_netcdf_var<std::string>(fo,"INSTRUME","TolTEC");
1146 add_netcdf_var(fo, "HWPR", calib.run_hwpr);
1147 add_netcdf_var<std::string>(fo, "TELESCOP", "LMT");
1148 add_netcdf_var<std::string>(fo, "PIPELINE", "CITLALI");
1149 add_netcdf_var<std::string>(fo, "VERSION", CITLALI_GIT_VERSION);
1150 add_netcdf_var<std::string>(fo, "KIDS", KIDSCPP_GIT_VERSION);
1151 add_netcdf_var<std::string>(fo, "TULA", TULA_GIT_VERSION);
1152 add_netcdf_var<std::string>(fo, "PROJID", telescope.project_id);
1153 add_netcdf_var<std::string>(fo, "GOAL", redu_type);
1154 add_netcdf_var<std::string>(fo, "OBSGOAL", telescope.obs_goal);
1155 add_netcdf_var<std::string>(fo, "TYPE", tod_type);
1156 add_netcdf_var<std::string>(fo, "GROUPING", map_grouping);
1157 add_netcdf_var<std::string>(fo, "METHOD", map_method);
1158 add_netcdf_var(fo, "EXPTIME", omb.exposure_time);
1159 add_netcdf_var<std::string>(fo, "RADESYS", telescope.pixel_axes);
1160 add_netcdf_var(fo, "TAN_RA", telescope.tel_header["Header.Source.Ra"][0]);
1161 add_netcdf_var(fo, "TAN_DEC", telescope.tel_header["Header.Source.Dec"][0]);
1162 add_netcdf_var(fo, "MEAN_EL", RAD_TO_DEG*telescope.tel_data["TelElAct"].mean());
1163 add_netcdf_var(fo, "MEAN_AZ", RAD_TO_DEG*telescope.tel_data["TelAzAct"].mean());
1164 add_netcdf_var(fo, "MEAN_PA", RAD_TO_DEG*telescope.tel_data["ActParAng"].mean());
1165
1166 // add beamsizes
1167 for (const auto &arr: calib.arrays) {
1168 if (std::get<0>(calib.array_fwhms[arr]) >= std::get<1>(calib.array_fwhms[arr])) {
1169 add_netcdf_var(fo, "BMAJ_"+toltec_io.array_name_map[arr], std::get<0>(calib.array_fwhms[arr]));
1170 add_netcdf_var(fo, "BMIN_"+toltec_io.array_name_map[arr], std::get<1>(calib.array_fwhms[arr]));
1172 }
1173 else {
1174 add_netcdf_var(fo, "BMAJ_"+toltec_io.array_name_map[arr], std::get<1>(calib.array_fwhms[arr]));
1175 add_netcdf_var(fo, "BMIN_"+toltec_io.array_name_map[arr], std::get<0>(calib.array_fwhms[arr]));
1176 add_netcdf_var(fo, "BPA_"+toltec_io.array_name_map[arr], (calib.array_pas[arr] + pi/2)*RAD_TO_DEG);
1177 }
1178 }
1179
1180 add_netcdf_var(fo, "BUNIT", omb.sig_unit);
1181
1182 // add jinc shape params
1183 if (map_method=="jinc") {
1184 add_netcdf_var(fo, "JINC_R", jinc_mm.r_max);
1185 for (const auto &[key,val]: toltec_io.array_name_map) {
1186 add_netcdf_var(fo, "JINC_A_"+val, jinc_mm.shape_params[calib.arrays(key)][0]);
1187 add_netcdf_var(fo, "JINC_B_"+val, jinc_mm.shape_params[calib.arrays(key)][0]);
1188 add_netcdf_var(fo, "JINC_C_"+val, jinc_mm.shape_params[calib.arrays(key)][0]);
1189 }
1190 }
1191
1192 // add mean tau
1193 if (rtcproc.run_extinction) {
1194 Eigen::VectorXd tau_el(1);
1195 tau_el << telescope.tel_data["TelElAct"].mean();
1196 auto tau_freq = rtcproc.calibration.calc_tau(tau_el, telescope.tau_225_GHz);
1197
1198 Eigen::Index i = 0;
1199 for (auto const& [key, val] : tau_freq) {
1200 add_netcdf_var(fo, "MEAN_TAU_"+toltec_io.array_name_map[calib.arrays(i)], val[0]);
1201 i++;
1202 }
1203 }
1204 else {
1205 for (Eigen::Index i=0; i<calib.arrays.size(); ++i) {
1206 add_netcdf_var(fo, "MEAN_TAU_"+toltec_io.array_name_map[calib.arrays(i)], 0.);
1207 }
1208 }
1209
1210 // add sample rate
1211 add_netcdf_var(fo, "SAMPRATE", telescope.fsmp);
1212
1213 // add apt table
1214 std::vector<string> apt_filename;
1215 std::stringstream ss(calib.apt_filepath);
1216 std::string item;
1217 char delim = '/';
1218
1219 while (getline (ss, item, delim)) {
1220 apt_filename.push_back(item);
1221 }
1222 add_netcdf_var<std::string>(fo, "APT", apt_filename.back());
1223
1224 add_netcdf_var(fo, "FRUITLOOPS_ITER", fruit_iter);
1225
1226 // add control/runtime parameters
1227 add_netcdf_var(fo, "CONFIG.VERBOSE", verbose_mode);
1228 add_netcdf_var(fo, "CONFIG.POLARIZED", rtcproc.run_polarization);
1229 add_netcdf_var(fo, "CONFIG.DESPIKED", rtcproc.run_despike);
1230 add_netcdf_var(fo, "CONFIG.TODFILTERED", rtcproc.run_tod_filter);
1231 add_netcdf_var(fo, "CONFIG.DOWNSAMPLED", rtcproc.run_downsample);
1232 add_netcdf_var(fo, "CONFIG.CALIBRATED", rtcproc.run_calibrate);
1233 add_netcdf_var(fo, "CONFIG.EXTINCTION", rtcproc.run_extinction);
1234 add_netcdf_var<std::string>(fo, "CONFIG.EXTINCTION.EXTMODEL", rtcproc.calibration.extinction_model);
1235 add_netcdf_var<std::string>(fo, "CONFIG.WEIGHT.TYPE", ptcproc.weighting_type);
1236 add_netcdf_var(fo, "CONFIG.INV_VAR.RTC.WTLOW", rtcproc.lower_inv_var_factor);
1237 add_netcdf_var(fo, "CONFIG.INV_VAR.RTC.WTHIGH", rtcproc.upper_inv_var_factor);
1238 add_netcdf_var(fo, "CONFIG.INV_VAR.PTC.WTLOW", ptcproc.lower_inv_var_factor);
1239 add_netcdf_var(fo, "CONFIG.INV_VAR.PTC.WTHIGH", ptcproc.upper_inv_var_factor);
1240 add_netcdf_var(fo, "CONFIG.WEIGHT.PTC.WTLOW", ptcproc.lower_weight_factor);
1241 add_netcdf_var(fo, "CONFIG.WEIGHT.PTC.WTHIGH", ptcproc.upper_weight_factor);
1242 add_netcdf_var(fo, "CONFIG.WEIGHT.MEDWTFACTOR", ptcproc.med_weight_factor);
1243 add_netcdf_var(fo, "CONFIG.CLEANED", ptcproc.run_clean);
1244
1245 // loop through arrays and add number of eigenvalues removed
1246 for (Eigen::Index i=0; i<calib.arrays.size(); ++i) {
1247 if (ptcproc.run_clean) {
1248 add_netcdf_var(fo, "CONFIG.CLEANED.NEIG_"+toltec_io.array_name_map[calib.arrays(i)],
1250 }
1251 else {
1252 add_netcdf_var(fo, "CONFIG.CLEANED.NEIG_"+toltec_io.array_name_map[calib.arrays(i)], 0);
1253 }
1254 }
1255
1256 // fruit loops parameters
1257 add_netcdf_var(fo, "CONFIG.FRUITLOOPS", ptcproc.run_fruit_loops);
1258 add_netcdf_var<std::string>(fo, "CONFIG.FRUITLOOPS.PATH", ptcproc.fruit_loops_path);
1259 add_netcdf_var(fo, "CONFIG.FRUITLOOPS.S2N", ptcproc.fruit_loops_sig2noise);
1260 for (Eigen::Index i=0; i<calib.arrays.size(); ++i) {
1262 add_netcdf_var(fo, "CONFIG.FRUITLOOPS.FLUX_"+toltec_io.array_name_map[calib.arrays(i)], ptcproc.fruit_loops_flux(calib.arrays(i)));
1263 }
1264 else {
1265 add_netcdf_var(fo, "CONFIG.FRUITLOOPS.FLUX_"+toltec_io.array_name_map[calib.arrays(i)], 0);
1266 }
1267 }
1268
1269 add_netcdf_var(fo, "CONFIG.FRUITLOOPS.MAXITER", ptcproc.fruit_loops_iters);
1270
1271 fo.close();
1272 }
1273}
1274
1275template <engine_utils::toltecIO::ProdType prod_t>
1277 // name for std map
1278 std::string name;
1279 // subdirectory name
1280 std::string dir_name = obsnum_dir_name + "raw/";
1281
1282 // if config subdirectory name is specified, add it
1283 if (tod_output_subdir_name != "null") {
1284 dir_name = dir_name + tod_output_subdir_name + "/";
1285 }
1286
1287 // rtc tod output filename setup
1288 if constexpr (prod_t == engine_utils::toltecIO::rtc_timestream) {
1293
1294 tod_filename["rtc"] = filename + ".nc";
1295 name = "rtc";
1296 }
1297
1298 // ptc tod output filename setup
1299 else if constexpr (prod_t == engine_utils::toltecIO::ptc_timestream) {
1304
1305 tod_filename["ptc"] = filename + ".nc";
1306 name = "ptc";
1307 }
1308
1309 // create netcdf file
1310 netCDF::NcFile fo(tod_filename[name], netCDF::NcFile::replace);
1311
1312 // add tod output type to file
1313 netCDF::NcDim n_tod_output_type_dim = fo.addDim("n_tod_output_type",1);
1314 netCDF::NcVar tod_output_type_var = fo.addVar("tod_output_type",netCDF::ncString, n_tod_output_type_dim);
1315 const std::vector<size_t> tod_output_type_index = {0};
1316
1317 if constexpr (prod_t == engine_utils::toltecIO::rtc_timestream) {
1318 std::string tod_output_type_name = "rtc";
1319 tod_output_type_var.putVar(tod_output_type_index,tod_output_type_name);
1320 }
1321 else if constexpr (prod_t == engine_utils::toltecIO::ptc_timestream) {
1322 std::string tod_output_type_name = "ptc";
1323 tod_output_type_var.putVar(tod_output_type_index,tod_output_type_name);
1324
1325 // number of eigenvalues
1326 netCDF::NcDim n_eigs_dim = fo.addDim("n_eigs",ptcproc.cleaner.n_calc);
1327 }
1328
1329 // add obsnum
1330 netCDF::NcVar obsnum_v = fo.addVar("obsnum",netCDF::ncInt);
1331 obsnum_v.putAtt("units","N/A");
1332 int obsnum_int = std::stoi(obsnum);
1333 obsnum_v.putVar(&obsnum_int);
1334
1335 // add source ra
1336 netCDF::NcVar source_ra_v = fo.addVar("SourceRa",netCDF::ncDouble);
1337 source_ra_v.putAtt("units","rad");
1338 source_ra_v.putVar(&telescope.tel_header["Header.Source.Ra"](0));
1339
1340 // add source dec
1341 netCDF::NcVar source_dec_v = fo.addVar("SourceDec",netCDF::ncDouble);
1342 source_dec_v.putAtt("units","rad");
1343 source_dec_v.putVar(&telescope.tel_header["Header.Source.Dec"](0));
1344
1345 netCDF::NcDim n_pts_dim = fo.addDim("n_pts");
1346 netCDF::NcDim n_raw_scan_indices_dim = fo.addDim("n_raw_scan_indices", telescope.scan_indices.rows());
1347 netCDF::NcDim n_scan_indices_dim = fo.addDim("n_scan_indices", 2);
1348 netCDF::NcDim n_scans_dim = fo.addDim("n_scans", telescope.scan_indices.cols());
1349
1350 netCDF::NcDim n_dets_dim = fo.addDim("n_dets", calib.n_dets);
1351
1352 std::vector<netCDF::NcDim> dims = {n_pts_dim, n_dets_dim};
1353 std::vector<netCDF::NcDim> raw_scans_dims = {n_scans_dim, n_raw_scan_indices_dim};
1354 std::vector<netCDF::NcDim> scans_dims = {n_scans_dim, n_scan_indices_dim};
1355
1356 // raw file scan indices
1357 netCDF::NcVar raw_scan_indices_v = fo.addVar("raw_scan_indices",netCDF::ncInt, raw_scans_dims);
1358 raw_scan_indices_v.putAtt("units","N/A");
1359 raw_scan_indices_v.putVar(telescope.scan_indices.data());
1360
1361 // scan indices for data
1362 netCDF::NcVar scan_indices_v = fo.addVar("scan_indices",netCDF::ncInt, scans_dims);
1363 scan_indices_v.putAtt("units","N/A");
1364
1365 // signal
1366 netCDF::NcVar signal_v = fo.addVar("signal",netCDF::ncDouble, dims);
1367 signal_v.putAtt("units",omb.sig_unit);
1368
1369 // chunk sizes
1370 std::vector<std::size_t> chunkSizes;
1371 // set chunk mode
1372 netCDF::NcVar::ChunkMode chunkMode = netCDF::NcVar::nc_CHUNKED;
1373
1374 // set chunking to mean scan size and n_dets
1375 chunkSizes.push_back(((telescope.scan_indices.row(3) - telescope.scan_indices.row(2)).array() + 1).mean());
1376 chunkSizes.push_back(calib.n_dets);
1377
1378 // set signal chunking
1379 signal_v.setChunking(chunkMode, chunkSizes);
1380
1381 // flags
1382 netCDF::NcVar flags_v = fo.addVar("flags",netCDF::ncDouble, dims);
1383 flags_v.putAtt("units","N/A");
1384 flags_v.setChunking(chunkMode, chunkSizes);
1385
1386 // kernel
1387 if (rtcproc.run_kernel) {
1388 netCDF::NcVar kernel_v = fo.addVar("kernel",netCDF::ncDouble, dims);
1389 kernel_v.putAtt("units","N/A");
1390 kernel_v.setChunking(chunkMode, chunkSizes);
1391 }
1392
1393 // detector lat
1394 netCDF::NcVar det_lat_v = fo.addVar("det_lat",netCDF::ncDouble, dims);
1395 det_lat_v.putAtt("units","rad");
1396 det_lat_v.setChunking(chunkMode, chunkSizes);
1397
1398 // detector lon
1399 netCDF::NcVar det_lon_v = fo.addVar("det_lon",netCDF::ncDouble, dims);
1400 det_lon_v.putAtt("units","rad");
1401 det_lon_v.setChunking(chunkMode, chunkSizes);
1402
1403 // calc absolute pointing if in radec frame
1404 if (telescope.pixel_axes == "radec") {
1405 // detector absolute ra
1406 netCDF::NcVar det_ra_v = fo.addVar("det_ra",netCDF::ncDouble, dims);
1407 det_ra_v.putAtt("units","rad");
1408 det_ra_v.setChunking(chunkMode, chunkSizes);
1409
1410 // detector absolute dec
1411 netCDF::NcVar det_dec_v = fo.addVar("det_dec",netCDF::ncDouble, dims);
1412 det_dec_v.putAtt("units","rad");
1413 det_dec_v.setChunking(chunkMode, chunkSizes);
1414 }
1415
1416 // add apt table
1417 for (auto const& x: calib.apt) {
1418 netCDF::NcVar apt_v = fo.addVar("apt_" + x.first,netCDF::ncDouble, n_dets_dim);
1419 apt_v.putAtt("units",calib.apt_header_units[x.first]);
1420 }
1421
1422 // add telescope parameters
1423 for (auto const& x: telescope.tel_data) {
1424 netCDF::NcVar tel_data_v = fo.addVar(x.first,netCDF::ncDouble, n_pts_dim);
1425 tel_data_v.putAtt("units","rad");
1426 tel_data_v.setChunking(chunkMode, chunkSizes);
1427 }
1428
1429 // add pointing offset parameters
1430 for (auto const& x: pointing_offsets_arcsec) {
1431 logger->info("pointing_offsets_arcsec.second {} {}",x.first, x.second);
1432 netCDF::NcVar offsets_v = fo.addVar("pointing_offset_"+x.first,netCDF::ncDouble, n_pts_dim);
1433 offsets_v.putAtt("units","arcsec");
1434 offsets_v.setChunking(chunkMode, chunkSizes);
1435 }
1436
1437 // add weights
1438 if constexpr (prod_t == engine_utils::toltecIO::ptc_timestream) {
1439 std::vector<netCDF::NcDim> weight_dims = {n_scans_dim, n_dets_dim};
1440 netCDF::NcVar weights_v = fo.addVar("weights",netCDF::ncDouble, weight_dims);
1441 weights_v.putAtt("units","("+omb.sig_unit+")^-2");
1442 }
1443
1444 // add hwpr
1446 if (calib.run_hwpr) {
1447 netCDF::NcVar hwpr_v = fo.addVar("hwpr",netCDF::ncDouble, n_pts_dim);
1448 hwpr_v.putAtt("units","rad");
1449 }
1450 }
1451
1452 // add tel header
1453 netCDF::NcDim tel_header_dim = fo.addDim("tel_header_n_pts", 1);
1454 for (const auto &[key,val]: telescope.tel_header) {
1455 netCDF::NcVar tel_header_v = fo.addVar(key,netCDF::ncDouble, tel_header_dim);
1456 tel_header_v.putVar(&val(0));
1457 }
1458
1459 fo.close();
1460}
1461
1462//template <TCDataKind tc_t>
1464 logger->info("reduction info");
1465 logger->info("obsnum: {}", obsnum);
1466 logger->info("map buffer rows: {}", omb.n_rows);
1467 logger->info("map buffer cols: {}", omb.n_cols);
1468 logger->info("number of maps: {}", omb.signal.size());
1469 logger->info("map units: {}", omb.sig_unit);
1470 logger->info("polarized reduction: {}", rtcproc.run_polarization);
1471
1472 // total size of all maps
1473 double mb_size_total = 0;
1474
1475 // make a rough estimate of memory usage for obs map buffer
1476 double omb_size = 8*omb.n_rows*omb.n_cols*(omb.signal.size() + omb.weight.size() +
1477 omb.kernel.size() + omb.coverage.size())/1e9;
1478
1479 logger->info("estimated size of map buffer {} GB", omb_size);
1480
1481 mb_size_total = mb_size_total + omb_size;
1482
1483 // print info if coadd is requested
1484 if (run_coadd) {
1485 logger->info("coadd map buffer rows: {}", cmb.n_rows);
1486 logger->info("coadd map buffer cols: {}", cmb.n_cols);
1487
1488 // make a rough estimate of memory usage for coadd map buffer
1489 double cmb_size = 8*cmb.n_rows*cmb.n_cols*(cmb.signal.size() + cmb.weight.size() +
1490 cmb.kernel.size() + cmb.coverage.size())/1e9;
1491
1492 logger->info("estimated size of coadd buffer {} GB", cmb_size);
1493
1494 mb_size_total = mb_size_total + cmb_size;
1495
1496 // output info if coadd noise maps are requested
1497 if (run_noise) {
1498 logger->info("coadd map buffer noise maps: {}", cmb.n_noise);
1499 // make a rough estimate of memory usage for coadd noise maps
1500 double nmb_size = 8*cmb.n_rows*cmb.n_cols*cmb.noise.size()*cmb.n_noise/1e9;
1501 logger->info("estimated size of noise buffer {} GB", nmb_size);
1502 mb_size_total = mb_size_total + nmb_size;
1503 }
1504 }
1505 else {
1506 // output info if obs noise maps are requested
1507 if (run_noise) {
1508 logger->info("observation map buffer noise maps: {}", omb.n_noise);
1509 // make a rough estimate of memory usage for obs noise maps
1510 double nmb_size = 8*omb.n_rows*omb.n_cols*omb.noise.size()*omb.n_noise/1e9;
1511 logger->info("estimated size of noise buffer {} GB", nmb_size);
1512 mb_size_total = mb_size_total + nmb_size;
1513 }
1514 }
1515
1516 logger->info("estimated size of all maps {} GB", mb_size_total);
1517 logger->info("number of scans: {}",telescope.scan_indices.cols());
1518
1519 // test getting memory usage for fun
1520 /*struct sysinfo memInfo;
1521 long long totalPhysMem = memInfo.totalram;
1522 totalPhysMem *= memInfo.mem_unit;
1523
1524 logger->info("total physical memory available {} GB", (totalPhysMem/1024)/1e7);*/
1525 logger->info("physical memory used {} GB", engine_utils::get_phys_memory()/1e7);
1526}
1527
1528template <TCDataKind tc_t>
1530
1531 logger->debug("writing summary files for chunk {}",in.index.data);
1532
1533 std::string filename = "chunk_summary_" + std::to_string(in.index.data);
1534
1535 // write summary log file
1536 std::ofstream f;
1537 f.open (obsnum_dir_name+"/logs/" + filename + ".log");
1538
1539 f << "Summary file for scan " << in.index.data << "\n";
1540 f << "-Citlali version: " << CITLALI_GIT_VERSION << "\n";
1541 f << "-Kidscpp version: " << KIDSCPP_GIT_VERSION << "\n";
1542 f << "-Time of time chunk creation: " + in.creation_time + "\n";
1543 f << "-Time of file writing: " << engine_utils::current_date_time() << "\n";
1544
1545 f << "-Reduction type: " << redu_type << "\n";
1546 f << "-TOD type: " << tod_type << "\n";
1547 f << "-TOD unit: " << omb.sig_unit << "\n";
1548 f << "-TOD chunk type: " << in.name << "\n";
1549
1550 f << "-Calibrated: " << in.status.calibrated << "\n";
1551 f << "-Extinction Corrected: " << in.status.extinction_corrected << "\n";
1552 f << "-Demodulated: " << in.status.demodulated << "\n";
1553 f << "-Kernel Generated: " << in.status.kernel_generated << "\n";
1554 f << "-Despiked: " << in.status.despiked << "\n";
1555 f << "-TOD filtered: " << in.status.tod_filtered << "\n";
1556 f << "-Downsampled: " << in.status.downsampled << "\n";
1557 f << "-Cleaned: " << in.status.cleaned << "\n";
1558
1559 f << "-Scan length: " << in.scans.data.rows() << "\n";
1560
1561 f << "-Number of detectors: " << in.scans.data.cols() << "\n";
1562 f << "-Number of detectors flagged in APT table: " << (calib.apt["flag"].array()!=0).count() << "\n";
1563 f << "-Number of detectors flagged below weight limit: " << in.n_dets_low <<"\n";
1564 f << "-Number of detectors flagged above weight limit: " << in.n_dets_high << "\n";
1565 Eigen::Index n_flagged = in.n_dets_low + in.n_dets_high + (calib.apt["flag"].array()!=0).count();
1566 f << "-Number of detectors flagged: " << n_flagged << " (" << 100*float(n_flagged)/float(in.scans.data.cols()) << "%)\n";
1567
1568 f << "-NaNs found: " << in.scans.data.array().isNaN().count() << "\n";
1569 f << "-Infs found: " << in.scans.data.array().isInf().count() << "\n";
1570 f << "-Data min: " << in.scans.data.minCoeff() << " " << omb.sig_unit << "\n";
1571 f << "-Data max: " << in.scans.data.maxCoeff() << " " << omb.sig_unit << "\n";
1572 f << "-Data mean: " << in.scans.data.mean() << " " << omb.sig_unit << "\n";
1573 f << "-Data median: " << tula::alg::median(in.scans.data) << " " << omb.sig_unit << "\n";
1574 f << "-Data stddev: " << engine_utils::calc_std_dev(in.scans.data) << " " << omb.sig_unit << "\n";
1575
1576 if (in.status.kernel_generated) {
1577 f << "-Kernel max: " << in.kernel.data.maxCoeff() << " " << omb.sig_unit << "\n";
1578 }
1579
1580 f.close();
1581}
1582
1583template <typename map_buffer_t>
1584void Engine::write_map_summary(map_buffer_t &mb) {
1585
1586 logger->debug("writing map summary files");
1587
1588 std::string filename = "map_summary";
1589 std::ofstream f;
1590 f.open (obsnum_dir_name+"/logs/" + filename + ".log");
1591
1592 f << "Summary file for maps\n";
1593 f << "-Citlali version: " << CITLALI_GIT_VERSION << "\n";
1594 f << "-Kidscpp version: " << KIDSCPP_GIT_VERSION << "\n";
1595 f << "-Time of file writing: " << engine_utils::current_date_time() << "\n";
1596
1597 f << "-Reduction type: " << redu_type << "\n";
1598 f << "-Map type: " << tod_type << "\n";
1599 f << "-Map grouping: " << map_grouping << "\n";
1600 f << "-Rows: " << mb.n_rows << "\n";
1601 f << "-Cols: " << mb.n_cols << "\n";
1602 f << "-Number of maps: " << n_maps << "\n";
1603 f << "-Signal map unit: " << mb.sig_unit << "\n";
1604 f << "-Weight map unit: " << "1/(" + mb.sig_unit + ")^2" << "\n";
1605 f << "-Kernel maps generated: " << !mb.kernel.empty() << "\n";
1606 f << "-Coverage maps generated: " << !mb.coverage.empty() << "\n";
1607 f << "-Noise maps generated: " << !mb.noise.empty() << "\n";
1608 f << "-Number of noise maps: " << mb.noise.size() << "\n";
1609
1610 // map to count nans for all maps
1611 std::map<std::string,int> n_nans;
1612 n_nans["signal"] = 0;
1613 n_nans["weight"] = 0;
1614 n_nans["kernel"] = 0;
1615 n_nans["coverage"] = 0;
1616 n_nans["noise"] = 0;
1617
1618 // maps to hold infs for all maps
1619 std::map<std::string,int> n_infs;
1620 n_infs["signal"] = 0;
1621 n_infs["weight"] = 0;
1622 n_nans["kernel"] = 0;
1623 n_infs["coverage"] = 0;
1624 n_infs["noise"] = 0;
1625
1626 // loop through maps and count up nans and infs
1627 for (Eigen::Index i=0; i<mb.signal.size(); ++i) {
1628 n_nans["signal"] = n_nans["signal"] + mb.signal[i].array().isNaN().count();
1629 n_nans["weight"] = n_nans["weight"] + mb.weight[i].array().isNaN().count();
1630
1631 // check kernel for nans if requested
1632 if (!mb.kernel.empty()) {
1633 n_nans["kernel"] = n_nans["kernel"] + mb.kernel[i].array().isNaN().count();
1634 }
1635 // check coverage map for nans if available
1636 if (!mb.coverage.empty()) {
1637 n_nans["coverage"] = n_nans["coverage"] + mb.coverage[i].array().isNaN().count();
1638 }
1639
1640 n_infs["signal"] = n_infs["signal"] + mb.signal[i].array().isInf().count();
1641 n_infs["weight"] = n_infs["weight"] + mb.weight[i].array().isInf().count();
1642
1643 // check kernel for infs if requested
1644 if (!mb.kernel.empty()) {
1645 n_infs["kernel"] = n_infs["kernel"] + mb.kernel[i].array().isInf().count();
1646 }
1647 // check coverage map for infs if available
1648 if (!mb.coverage.empty()) {
1649 n_infs["coverage"] = n_infs["coverage"] + mb.coverage[i].array().isInf().count();
1650 }
1651
1652 // loop through noise maps and check for nans and infs
1653 if (!mb.noise.empty()) {
1654 for (Eigen::Index j=0; j<mb.noise.size(); ++j) {
1655 Eigen::Map<Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic>> noise_matrix(mb.noise[i].data() + j * mb.n_rows * mb.n_cols,
1656 mb.n_rows, mb.n_cols);
1657
1658 n_nans["noise"] = n_nans["noise"] + noise_matrix.array().isNaN().count();
1659 n_infs["noise"] = n_infs["noise"] + noise_matrix.array().isInf().count();
1660 }
1661 }
1662 }
1663
1664 for (auto const& [key, val] : n_nans) {
1665 f << "-Number of "+ key + " NaNs: " << val << "\n";
1666 }
1667
1668 for (auto const& [key, val] : n_infs) {
1669 f << "-Number of "+ key + " Infs: " << val << "\n";
1670 }
1671}
1672
1673template <mapmaking::MapType map_t, engine_utils::toltecIO::DataType data_t, engine_utils::toltecIO::ProdType prod_t>
1674auto Engine::setup_filenames(std::string dir_name) {
1675
1676 std::string filename;
1677
1678 // raw obs maps
1679 if constexpr (map_t == mapmaking::RawObs) {
1680 filename = toltec_io.create_filename<data_t, prod_t, engine_utils::toltecIO::raw>
1681 (dir_name, redu_type, "", obsnum, telescope.sim_obs);
1682 }
1683 // filtered obs maps
1684 else if constexpr (map_t == mapmaking::FilteredObs) {
1686 (dir_name, redu_type, "", obsnum, telescope.sim_obs);
1687 }
1688 // raw coadded maps
1689 else if constexpr (map_t == mapmaking::RawCoadd) {
1690 filename = toltec_io.create_filename<data_t, prod_t, engine_utils::toltecIO::raw>
1691 (dir_name, "", "", "", telescope.sim_obs);
1692 }
1693 // filtered coadded maps
1694 else if constexpr (map_t == mapmaking::FilteredCoadd) {
1696 (dir_name, "", "", "", telescope.sim_obs);
1697 }
1698
1699 return filename;
1700}
1701
1703 // get name for extension layer
1704 std::string map_name = "";
1705
1706 // only update name if we're not in array mode
1707 if (map_grouping!="array") {
1708 // if in nw mode
1709 if (map_grouping=="nw") {
1710 map_name = map_name + "nw_" + std::to_string(calib.nws(i)) + "_";
1711 }
1712 else if (map_grouping=="fg") {
1713 // find all detectors belonging to each fg
1714 Eigen::VectorXI array_indices(calib.fg.size()*calib.n_arrays*rtcproc.polarization.stokes_params.size());
1715 Eigen::Index k = 0;
1716 for (Eigen::Index j=0; j<calib.n_arrays; ++j) {
1717 for (Eigen::Index l=0; l<rtcproc.polarization.stokes_params.size(); ++l) {
1718 for (Eigen::Index m=0; m<calib.fg.size(); ++m) {
1719 array_indices(k) = calib.fg(m);
1720 k++;
1721 }
1722 }
1723 }
1724 // if in fg mode
1725 map_name = map_name + "fg_" + std::to_string(array_indices(i)) + "_";
1726 }
1727 // if in detector mode
1728 else if (map_grouping=="detector") {
1729 map_name = map_name + "det_" + std::to_string(i) + "_";
1730 }
1731 }
1732
1733 return map_name;
1734}
1735
1736template <typename fits_io_type, class map_buffer_t>
1737void Engine::add_phdu(fits_io_type &fits_io, map_buffer_t &mb, Eigen::Index i) {
1738 // array name
1739 std::string name = toltec_io.array_name_map[calib.arrays(i)];
1740
1741 logger->debug("adding unit conversions");
1742
1743 // conversion to uK
1744 auto fwhm = (std::get<0>(calib.array_fwhms[calib.arrays(i)]) + std::get<1>(calib.array_fwhms[calib.arrays(i)]))/2;
1746
1747 // beam area in steradians
1748 auto beam_area_rad = 2.*pi*pow(fwhm*FWHM_TO_STD*ASEC_TO_RAD,2);
1749 // get Jy/pixel
1750 auto mJy_beam_to_Jy_px = 1e-3/beam_area_rad*pow(mb->pixel_size_rad,2);
1751
1752 // add unit conversions
1753 if (rtcproc.run_calibrate) {
1754 if (mb->sig_unit == "mJy/beam") {
1755 // conversion to mJy/beam
1756 fits_io->at(i).pfits->pHDU().addKey("to_mJy/beam", 1, "Conversion to mJy/beam");
1757 // conversion to MJy/sr
1758 fits_io->at(i).pfits->pHDU().addKey("to_MJy/sr", 1/(calib.array_beam_areas[calib.arrays(i)]*MJY_SR_TO_mJY_ASEC),
1759 "Conversion to MJy/sr");
1760 // conversion to uK
1761 fits_io->at(i).pfits->pHDU().addKey("to_uK", mJy_beam_to_uK, "Conversion to uK");
1762 // conversion to Jy/pixel
1763 fits_io->at(i).pfits->pHDU().addKey("to_Jy/pixel", mJy_beam_to_Jy_px, "Conversion to Jy/pixel");
1764 }
1765 else if (mb->sig_unit == "MJy/sr") {
1766 // conversion to mJy/beam
1767 fits_io->at(i).pfits->pHDU().addKey("to_mJy/beam", calib.array_beam_areas[calib.arrays(i)]*MJY_SR_TO_mJY_ASEC,
1768 "Conversion to mJy/beam");
1769 // conversion to MJy/Sr
1770 fits_io->at(i).pfits->pHDU().addKey("to_MJy/sr", 1, "Conversion to MJy/sr");
1771 // conversion to uK
1772 fits_io->at(i).pfits->pHDU().addKey("to_uK", calib.array_beam_areas[calib.arrays(i)]*MJY_SR_TO_mJY_ASEC*mJy_beam_to_uK,
1773 "Conversion to uK");
1774 // conversion to Jy/pixel
1775 fits_io->at(i).pfits->pHDU().addKey("to_Jy/pixel", calib.array_beam_areas[calib.arrays(i)]*MJY_SR_TO_mJY_ASEC*mJy_beam_to_Jy_px,
1776 "Conversion to Jy/pixel");
1777 }
1778 else if (mb->sig_unit == "uK") {
1779 // conversion to mJy/beam
1780 fits_io->at(i).pfits->pHDU().addKey("to_mJy/beam", 1/mJy_beam_to_uK, "Conversion to mJy/beam");
1781 // conversion to MJy/sr
1782 fits_io->at(i).pfits->pHDU().addKey("to_MJy/sr", 1/mJy_beam_to_uK/(calib.array_beam_areas[calib.arrays(i)]*MJY_SR_TO_mJY_ASEC),
1783 "Conversion to MJy/sr");
1784 // conversion to uK
1785 fits_io->at(i).pfits->pHDU().addKey("to_uK", 1, "Conversion to uK");
1786 // conversion to Jy/pixel
1787 fits_io->at(i).pfits->pHDU().addKey("to_Jy/pixel", (1/mJy_beam_to_uK)*mJy_beam_to_Jy_px, "Conversion to Jy/pixel");
1788 }
1789 else if (mb->sig_unit == "Jy/pixel") {
1790 // conversion to mJy/beam
1791 fits_io->at(i).pfits->pHDU().addKey("to_mJy/beam", 1/mJy_beam_to_Jy_px, "Conversion to mJy/beam");
1792 // conversion to MJy/sr
1793 fits_io->at(i).pfits->pHDU().addKey("to_MJy/sr", (1/mJy_beam_to_Jy_px)/(calib.array_beam_areas[calib.arrays(i)]*MJY_SR_TO_mJY_ASEC),
1794 "Conversion to MJy/sr");
1795 // conversion to uK
1796 fits_io->at(i).pfits->pHDU().addKey("to_uK", mJy_beam_to_uK/mJy_beam_to_Jy_px, "Conversion to uK");
1797 // conversion to Jy/pixel
1798 fits_io->at(i).pfits->pHDU().addKey("to_Jy/pixel", 1, "Conversion to Jy/pixel");
1799 }
1800 }
1801 // if flux calibration is disabled
1802 else {
1803 fits_io->at(i).pfits->pHDU().addKey("to_mJy/beam", "N/A", "Conversion to mJy/beam");
1804 fits_io->at(i).pfits->pHDU().addKey("to_MJy/sr", "N/A", "Conversion to MJy/sr");
1805 fits_io->at(i).pfits->pHDU().addKey("to_uK", "N/A", "Conversion to uK");
1806 fits_io->at(i).pfits->pHDU().addKey("to_Jy/pixel", "N/A", "Conversion to Jy/pixel");
1807 }
1808
1809 // add source flux for beammaps
1810 if (redu_type == "beammap") {
1811 fits_io->at(i).pfits->pHDU().addKey("HEADER.SOURCE.FLUX_MJYPERBEAM", beammap_fluxes_mJy_beam[name], "Source flux (mJy/beam)");
1812 fits_io->at(i).pfits->pHDU().addKey("HEADER.SOURCE.FLUX_MJYPERSR", beammap_fluxes_MJy_Sr[name], "Source flux (MJy/sr)");
1813
1814 fits_io->at(i).pfits->pHDU().addKey("BEAMMAP.ITER_TOLERANCE", beammap_iter_tolerance, "Beammap iteration tolerance");
1815 fits_io->at(i).pfits->pHDU().addKey("BEAMMAP.ITER_MAX", beammap_iter_max, "Beammap max iterations");
1816 fits_io->at(i).pfits->pHDU().addKey("BEAMMAP.IS_DEROTATED", beammap_derotate, "Beammap derotated");
1817 // add reference detector information
1819 fits_io->at(i).pfits->pHDU().addKey("BEAMMAP.REF_DET_INDEX", beammap_reference_det, "Beammap Reference det (rotation center)");
1820 fits_io->at(i).pfits->pHDU().addKey("BEAMMAP.REF_X_T", calib.apt["x_t"](beammap_reference_det), "Az rotation center (arcsec)");
1821 fits_io->at(i).pfits->pHDU().addKey("BEAMMAP.REF_Y_T", calib.apt["y_t"](beammap_reference_det), "Alt rotation center (arcsec)");
1822 }
1823 else {
1824 fits_io->at(i).pfits->pHDU().addKey("BEAMMAP.REF_DET_INDEX", -99, "Beammap Reference det (rotation center)");
1825 fits_io->at(i).pfits->pHDU().addKey("BEAMMAP.REF_X_T", "N/A", "Az rotation center (arcsec)");
1826 fits_io->at(i).pfits->pHDU().addKey("BEAMMAP.REF_Y_T", "N/A", "Alt rotation center (arcsec)");
1827 }
1828 }
1829
1830 logger->debug("adding obsnums");
1831
1832 // add obsnums
1833 for (Eigen::Index j=0; j<mb->obsnums.size(); ++j) {
1834 fits_io->at(i).pfits->pHDU().addKey("OBSNUM"+std::to_string(j), mb->obsnums.at(j), "Observation Number " + std::to_string(j));
1835 }
1836
1837 // add date and time of obs
1838 if (mb->obsnums.size()==1) {
1839 fits_io->at(i).pfits->pHDU().addKey("DATEOBS0", date_obs.back(), "Date and time of observation 0");
1840 }
1841 else {
1842 for (Eigen::Index j=0; j<mb->obsnums.size(); ++j) {
1843 fits_io->at(i).pfits->pHDU().addKey("DATEOBS"+std::to_string(j), date_obs[j], "Date and time of observation "+std::to_string(j));
1844 }
1845 }
1846
1847 logger->debug("adding obs info");
1848
1849 // add source
1850 fits_io->at(i).pfits->pHDU().addKey("SOURCE", telescope.source_name, "Source name");
1851 // add instrument
1852 fits_io->at(i).pfits->pHDU().addKey("INSTRUME", "TolTEC", "Instrument");
1853 // add hwpr
1854 fits_io->at(i).pfits->pHDU().addKey("HWPR", calib.run_hwpr, "HWPR installed");
1855 // add telescope
1856 fits_io->at(i).pfits->pHDU().addKey("TELESCOP", "LMT", "Telescope");
1857 // add wavelength
1858 fits_io->at(i).pfits->pHDU().addKey("WAV", name, "Wavelength");
1859 // add pipeline
1860 fits_io->at(i).pfits->pHDU().addKey("PIPELINE", "CITLALI", "Redu pipeline");
1861 // add citlali version
1862 fits_io->at(i).pfits->pHDU().addKey("VERSION", CITLALI_GIT_VERSION, "CITLALI_GIT_VERSION");
1863 // add kids version
1864 fits_io->at(i).pfits->pHDU().addKey("KIDS", KIDSCPP_GIT_VERSION, "KIDSCPP_GIT_VERSION");
1865 // add kids version
1866 fits_io->at(i).pfits->pHDU().addKey("TULA", TULA_GIT_VERSION, "TULA_GIT_VERSION");
1867 // project id
1868 fits_io->at(i).pfits->pHDU().addKey("PROJID", telescope.project_id, "Project ID");
1869 // add redu type
1870 fits_io->at(i).pfits->pHDU().addKey("GOAL", redu_type, "Reduction type");
1871 // add obs goal
1872 fits_io->at(i).pfits->pHDU().addKey("OBSGOAL", telescope.obs_goal, "Obs goal");
1873 // add tod type
1874 fits_io->at(i).pfits->pHDU().addKey("TYPE", tod_type, "TOD Type");
1875 // add map grouping
1876 fits_io->at(i).pfits->pHDU().addKey("GROUPING", map_grouping, "Map grouping");
1877 // add map grouping
1878 fits_io->at(i).pfits->pHDU().addKey("METHOD", map_method, "Map method");
1879 // add exposure time
1880 fits_io->at(i).pfits->pHDU().addKey("EXPTIME", mb->exposure_time, "Exposure time (sec)");
1881 // add pixel axes
1882 fits_io->at(i).pfits->pHDU().addKey("RADESYS", telescope.pixel_axes, "Coord Reference Frame");
1883 // add source ra
1884 fits_io->at(i).pfits->pHDU().addKey("SRC_RA", telescope.tel_header["Header.Source.Ra"][0], "Source RA (radians)");
1885 // add source dec
1886 fits_io->at(i).pfits->pHDU().addKey("SRC_DEC", telescope.tel_header["Header.Source.Dec"][0], "Source Dec (radians)");
1887 // add map tangent point ra
1888 fits_io->at(i).pfits->pHDU().addKey("TAN_RA", telescope.tel_header["Header.Source.Ra"][0], "Map Tangent Point RA (radians)");
1889 //add map tangent point dec
1890 fits_io->at(i).pfits->pHDU().addKey("TAN_DEC", telescope.tel_header["Header.Source.Dec"][0], "Map Tangent Point Dec (radians)");
1891 // add mean alt
1892 fits_io->at(i).pfits->pHDU().addKey("MEAN_EL", RAD_TO_DEG*telescope.tel_data["TelElAct"].mean(), "Mean Elevation (deg)");
1893 // add mean az
1894 fits_io->at(i).pfits->pHDU().addKey("MEAN_AZ", RAD_TO_DEG*telescope.tel_data["TelAzAct"].mean(), "Mean Azimuth (deg)");
1895 // add mean parallactic angle
1896 fits_io->at(i).pfits->pHDU().addKey("MEAN_PA", RAD_TO_DEG*telescope.tel_data["ActParAng"].mean(), "Mean Parallactic angle (deg)");
1897
1898 logger->debug("adding beamsizes");
1899
1900 // add beamsizes
1901 if (std::get<0>(calib.array_fwhms[calib.arrays(i)]) >= std::get<1>(calib.array_fwhms[calib.arrays(i)])) {
1902 fits_io->at(i).pfits->pHDU().addKey("BMAJ", std::get<0>(calib.array_fwhms[calib.arrays(i)]), "beammaj (arcsec)");
1903 fits_io->at(i).pfits->pHDU().addKey("BMIN", std::get<1>(calib.array_fwhms[calib.arrays(i)]), "beammin (arcsec)");
1904 fits_io->at(i).pfits->pHDU().addKey("BPA", calib.array_pas[calib.arrays(i)]*RAD_TO_DEG, "beampa (deg)");
1905 }
1906 else {
1907 fits_io->at(i).pfits->pHDU().addKey("BMAJ", std::get<1>(calib.array_fwhms[calib.arrays(i)]), "beammaj (arcsec)");
1908 fits_io->at(i).pfits->pHDU().addKey("BMIN", std::get<0>(calib.array_fwhms[calib.arrays(i)]), "beammin (arcsec)");
1909 fits_io->at(i).pfits->pHDU().addKey("BPA", (calib.array_pas[calib.arrays(i)] + pi/2)*RAD_TO_DEG, "beampa (deg)");
1910 }
1911
1912 fits_io->at(i).pfits->pHDU().addKey("BUNIT", mb->sig_unit, "bunit");
1913
1914 // add jinc shape params
1915 if (map_method=="jinc") {
1916 logger->debug("adding jinc params");
1917
1918 fits_io->at(i).pfits->pHDU().addKey("JINC_R", jinc_mm.r_max, "Jinc filter R_max");
1919 fits_io->at(i).pfits->pHDU().addKey("JINC_A", jinc_mm.shape_params[calib.arrays(i)][0], "Jinc filter param a");
1920 fits_io->at(i).pfits->pHDU().addKey("JINC_B", jinc_mm.shape_params[calib.arrays(i)][1], "Jinc filter param b");
1921 fits_io->at(i).pfits->pHDU().addKey("JINC_C", jinc_mm.shape_params[calib.arrays(i)][2], "Jinc filter param c");
1922 }
1923
1924 // add mean tau
1925 logger->debug("adding extinction");
1926 if (rtcproc.run_extinction) {
1927 Eigen::VectorXd tau_el(1);
1928 tau_el << telescope.tel_data["TelElAct"].mean();
1929 auto tau_freq = rtcproc.calibration.calc_tau(tau_el, telescope.tau_225_GHz);
1930
1931 fits_io->at(i).pfits->pHDU().addKey("MEAN_TAU", tau_freq[i](0), "mean tau (" + name + ")");
1932 }
1933 else {
1934 fits_io->at(i).pfits->pHDU().addKey("MEAN_TAU", 0., "mean tau (" + name + ")");
1935 }
1936
1937 // add sample rate
1938 fits_io->at(i).pfits->pHDU().addKey("SAMPRATE", telescope.fsmp, "sample rate (Hz)");
1939
1940 // add apt table to header
1941 if (mb->obsnums.size()==1) {
1942 std::vector<string> apt_filename;
1943 std::stringstream ss(calib.apt_filepath);
1944 std::string item;
1945 char delim = '/';
1946
1947 while (getline (ss, item, delim)) {
1948 apt_filename.push_back(item);
1949 }
1950 fits_io->at(i).pfits->pHDU().addKey("APT", apt_filename.back(), "APT table used");
1951 }
1952
1953 double rms;
1954
1955 if (redu_type != "beammap") {
1956 // estimate rms from weight maps
1957 rms = pow(mb->median_err(i),0.5);
1958 }
1959 else {
1960 rms = 0.0;
1961 }
1962
1963 // out-of-focus holography parameters
1964 if (! telescope.sim_obs) {
1965 logger->debug("adding oof params");
1966 fits_io->at(i).pfits->pHDU().addKey("OOF_RMS", rms, "rms of map background (" + mb->sig_unit +")");
1967 fits_io->at(i).pfits->pHDU().addKey("OOF_W", toltec_io.array_wavelength_map[calib.arrays(i)]/1000., "wavelength (m)");
1968 fits_io->at(i).pfits->pHDU().addKey("OOF_ID", static_cast<int>(toltec_io.array_wavelength_map[calib.arrays(i)]*1000), "instrument id");
1969 fits_io->at(i).pfits->pHDU().addKey("OOF_T", 3.0, "taper (dB)");
1970 fits_io->at(i).pfits->pHDU().addKey("OOF_M2X", telescope.tel_header["Header.M2.XReq"](0)/1000.*1e6, "oof m2x (microns)");
1971 fits_io->at(i).pfits->pHDU().addKey("OOF_M2Y", telescope.tel_header["Header.M2.YReq"](0)/1000.*1e6, "oof m2y (microns)");
1972 fits_io->at(i).pfits->pHDU().addKey("OOF_M2Z", telescope.tel_header["Header.M2.ZReq"](0)/1000.*1e6, "oof m2z (microns)");
1973
1974 fits_io->at(i).pfits->pHDU().addKey("OOF_RO", 25., "outer diameter of the antenna (m)");
1975 fits_io->at(i).pfits->pHDU().addKey("OOF_RI", 1.65, "inner diameter of the antenna (m)");
1976 }
1977
1978 fits_io->at(i).pfits->pHDU().addKey("FRUITLOOPS_ITER", fruit_iter, "Current fruit loops iteration");
1979
1980 // add control/runtime parameters
1981 logger->debug("adding config params");
1982 fits_io->at(i).pfits->pHDU().addKey("CONFIG.VERBOSE", verbose_mode, "Reduced in verbose mode");
1983 fits_io->at(i).pfits->pHDU().addKey("CONFIG.POLARIZED", rtcproc.run_polarization, "Polarized Obs");
1984 fits_io->at(i).pfits->pHDU().addKey("CONFIG.DESPIKED", rtcproc.run_despike, "Despiked");
1985 fits_io->at(i).pfits->pHDU().addKey("CONFIG.TODFILTERED", rtcproc.run_tod_filter, "TOD Filtered");
1986 fits_io->at(i).pfits->pHDU().addKey("CONFIG.DOWNSAMPLED", rtcproc.run_downsample, "Downsampled");
1987 fits_io->at(i).pfits->pHDU().addKey("CONFIG.CALIBRATED", rtcproc.run_calibrate, "Calibrated");
1988 fits_io->at(i).pfits->pHDU().addKey("CONFIG.EXTINCTION", rtcproc.run_extinction, "Extinction corrected");
1989 fits_io->at(i).pfits->pHDU().addKey("CONFIG.EXTINCTION.EXTMODEL", rtcproc.calibration.extinction_model, "Extinction model");
1990 fits_io->at(i).pfits->pHDU().addKey("CONFIG.WEIGHT.TYPE", ptcproc.weighting_type, "Weighting scheme");
1991 fits_io->at(i).pfits->pHDU().addKey("CONFIG.INV_VAR.RTC.WTLOW", rtcproc.lower_inv_var_factor, "RTC lower inv var cutoff");
1992 fits_io->at(i).pfits->pHDU().addKey("CONFIG.INV_VAR.RTC.WTHIGH", rtcproc.upper_inv_var_factor, "RTC upper inv var cutoff");
1993 fits_io->at(i).pfits->pHDU().addKey("CONFIG.INV_VAR.PTC.WTLOW", ptcproc.lower_inv_var_factor, "PTC lower inv var cutoff");
1994 fits_io->at(i).pfits->pHDU().addKey("CONFIG.INV_VAR.PTC.WTHIGH", ptcproc.upper_inv_var_factor, "PTC upper inv var cutoff");
1995 fits_io->at(i).pfits->pHDU().addKey("CONFIG.WEIGHT.PTC.WTLOW", ptcproc.lower_weight_factor, "PTC lower weight cutoff");
1996 fits_io->at(i).pfits->pHDU().addKey("CONFIG.WEIGHT.PTC.WTHIGH", ptcproc.upper_weight_factor, "PTC upper weight cutoff");
1997 fits_io->at(i).pfits->pHDU().addKey("CONFIG.WEIGHT.MEDWTFACTOR", ptcproc.med_weight_factor, "Median weight factor");
1998 fits_io->at(i).pfits->pHDU().addKey("CONFIG.CLEANED", ptcproc.run_clean, "Cleaned");
1999 if (ptcproc.run_clean) {
2000 fits_io->at(i).pfits->pHDU().addKey("CONFIG.CLEANED.NEIG", ptcproc.cleaner.n_eig_to_cut[calib.arrays(i)].sum(),
2001 "Number of eigenvalues removed");
2002 }
2003 else {
2004 fits_io->at(i).pfits->pHDU().addKey("CONFIG.CLEANED.NEIG", 0, "Number of eigenvalues removed");
2005 }
2006
2007 fits_io->at(i).pfits->pHDU().addKey("CONFIG.FRUITLOOPS", ptcproc.run_fruit_loops, "Fruit loops");
2008 fits_io->at(i).pfits->pHDU().addKey("CONFIG.FRUITLOOPS.PATH", ptcproc.fruit_loops_path, "Fruit loops path");
2009 fits_io->at(i).pfits->pHDU().addKey("CONFIG.FRUITLOOPS.TYPE", ptcproc.fruit_loops_type, "Fruit loops type");
2010 fits_io->at(i).pfits->pHDU().addKey("CONFIG.FRUITLOOPS.S2N", ptcproc.fruit_loops_sig2noise, "Fruit loops S/N");
2012 fits_io->at(i).pfits->pHDU().addKey("CONFIG.FRUITLOOPS.FLUX", ptcproc.fruit_loops_flux(calib.arrays(i)), "Fruit loops flux (" + mb->sig_unit + ")");
2013 }
2014 else {
2015 fits_io->at(i).pfits->pHDU().addKey("CONFIG.FRUITLOOPS.FLUX", 0, "Fruit loops flux (" + mb->sig_unit + ")");
2016 }
2017 fits_io->at(i).pfits->pHDU().addKey("CONFIG.FRUITLOOPS.MAXITER", ptcproc.fruit_loops_iters, "Fruit loops iterations");
2018
2019 // add telescope file header information
2020 if (mb->obsnums.size()==1) {
2021 logger->debug("adding tel params");
2022 for (auto const& [key, val] : telescope.tel_header) {
2023 logger->debug("adding {}: {}", key, val);
2024 fits_io->at(i).pfits->pHDU().addKey(key, val(0), key);
2025 }
2026 }
2027}
2028
2029template <typename fits_io_type, class map_buffer_t>
2030void Engine::write_maps(fits_io_type &fits_io, fits_io_type &noise_fits_io, map_buffer_t &mb, Eigen::Index i) {
2031 // get name for extension layer
2032 std::string map_name = get_map_name(i);
2033
2034 // get the array for the given map
2035 Eigen::Index map_index = arrays_to_maps(i);
2036 // get the stokes parameter for the given map
2037 Eigen::Index stokes_index = maps_to_stokes(i);
2038
2039 // update wcs ctypes for frequency and stokes params
2040 mb->wcs.crval[2] = toltec_io.array_freq_map[calib.arrays[maps_to_arrays(i)]];
2041 mb->wcs.crval[3] = stokes_index;
2042
2043 // signal map
2044 fits_io->at(map_index).add_hdu("signal_" + map_name + rtcproc.polarization.stokes_params[stokes_index], mb->signal[i]);
2045 fits_io->at(map_index).add_wcs(fits_io->at(map_index).hdus.back(), mb->wcs, telescope.tel_header["Header.Source.Epoch"](0));
2046 fits_io->at(map_index).hdus.back()->addKey("UNIT", mb->sig_unit, "Unit of map");
2047
2048 // weight map
2049 fits_io->at(map_index).add_hdu("weight_" + map_name + rtcproc.polarization.stokes_params[stokes_index], mb->weight[i]);
2050 fits_io->at(map_index).add_wcs(fits_io->at(map_index).hdus.back(), mb->wcs, telescope.tel_header["Header.Source.Epoch"](0));
2051 fits_io->at(map_index).hdus.back()->addKey("UNIT", "1/("+mb->sig_unit+")^2", "Unit of map");
2052 if (redu_type != "beammap" && std::fabs(mb->median_err(i)) > std::numeric_limits<double>::epsilon()) {
2053 fits_io->at(map_index).hdus.back()->addKey("MEDERR", pow(mb->median_err(i),0.5), "Median Error ("+mb->sig_unit+")");
2054 }
2055 else {
2056 fits_io->at(map_index).hdus.back()->addKey("MEDERR", 0.0, "Median Error ("+mb->sig_unit+")");
2057 }
2058
2059 // kernel map
2060 if (rtcproc.run_kernel) {
2061 fits_io->at(map_index).add_hdu("kernel_" + map_name + rtcproc.polarization.stokes_params[stokes_index], mb->kernel[i]);
2062 fits_io->at(map_index).hdus.back()->addKey("TYPE",rtcproc.kernel.type, "Kernel type");
2063
2064 // add fwhm
2065 double fwhm = -99;
2066 if (rtcproc.kernel.type!="fits") {
2067 if (rtcproc.kernel.fwhm_rad<=0) {
2068 fwhm = (std::get<0>(calib.array_fwhms[calib.arrays(i)]) + std::get<1>(calib.array_fwhms[calib.arrays(i)]))/2;
2069 }
2070 else {
2072 }
2073 }
2074 fits_io->at(map_index).hdus.back()->addKey("FWHM",fwhm,"Kernel fwhm (arcsec)");
2075 fits_io->at(map_index).add_wcs(fits_io->at(map_index).hdus.back(), mb->wcs, telescope.tel_header["Header.Source.Epoch"](0));
2076 fits_io->at(map_index).hdus.back()->addKey("UNIT", mb->sig_unit, "Unit of map");
2077 }
2078
2079 // coverage map
2080 if (!mb->coverage.empty()) {
2081 fits_io->at(map_index).add_hdu("coverage_" + map_name + rtcproc.polarization.stokes_params[stokes_index], mb->coverage[i]);
2082 fits_io->at(map_index).add_wcs(fits_io->at(map_index).hdus.back(), mb->wcs, telescope.tel_header["Header.Source.Epoch"](0));
2083 fits_io->at(map_index).hdus.back()->addKey("UNIT", "sec", "Unit of map");
2084 }
2085
2086 /* coverage bool and signal-to-noise maps */
2087 if (!mb->coverage.empty()) {
2088 // need these to use eigen select
2089 Eigen::MatrixXd ones, zeros;
2090 ones.setOnes(mb->weight[i].rows(), mb->weight[i].cols());
2091 zeros.setZero(mb->weight[i].rows(), mb->weight[i].cols());
2092
2093 // get weight threshold for current map
2094 auto [weight_threshold, cov_ranges, cov_n_rows, cov_n_cols] = mb->calc_cov_region(i);
2095 // if weight is less than threshold, set to zero, otherwise set to one
2096 Eigen::MatrixXd coverage_bool = (mb->weight[i].array() < weight_threshold).select(zeros,ones);
2097
2098 // coverage bool map
2099 fits_io->at(map_index).add_hdu("coverage_bool_" + map_name + rtcproc.polarization.stokes_params[stokes_index], coverage_bool);
2100 fits_io->at(map_index).add_wcs(fits_io->at(map_index).hdus.back(), mb->wcs, telescope.tel_header["Header.Source.Epoch"](0));
2101 fits_io->at(map_index).hdus.back()->addKey("UNIT", "N/A", "Unit of map");
2102 fits_io->at(map_index).hdus.back()->addKey("WTTHRESH", weight_threshold, "Weight threshold");
2103
2104 // signal-to-noise map
2105 Eigen::MatrixXd sig2noise = mb->signal[i].array()*sqrt(mb->weight[i].array());
2106 fits_io->at(map_index).add_hdu("sig2noise_" + map_name + rtcproc.polarization.stokes_params[stokes_index], sig2noise);
2107 fits_io->at(map_index).add_wcs(fits_io->at(map_index).hdus.back(), mb->wcs, telescope.tel_header["Header.Source.Epoch"](0));
2108 fits_io->at(map_index).hdus.back()->addKey("UNIT", "N/A", "Unit of map");
2109 }
2110
2111 // write noise maps
2112 if (!mb->noise.empty()) {
2113 for (Eigen::Index n=0; n<mb->n_noise; ++n) {
2114 Eigen::Map<Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic>> noise_matrix(mb->noise[i].data() + n * mb->n_rows * mb->n_cols,
2115 mb->n_rows, mb->n_cols);
2116
2117 noise_fits_io->at(map_index).add_hdu("signal_" + map_name + std::to_string(n) + "_" + rtcproc.polarization.stokes_params[stokes_index],
2118 noise_matrix);
2119 noise_fits_io->at(map_index).add_wcs(noise_fits_io->at(map_index).hdus.back(), mb->wcs, telescope.tel_header["Header.Source.Epoch"](0));
2120 noise_fits_io->at(map_index).hdus.back()->addKey("UNIT", mb->sig_unit, "Unit of map");
2121 noise_fits_io->at(map_index).hdus.back()->addKey("MEDRMS", mb->median_rms[i], "Median RMS of noise maps");
2122 }
2123 }
2124}
2125
2126template <mapmaking::MapType map_t, class map_buffer_t>
2127void Engine::write_psd(map_buffer_t &mb, std::string dir_name) {
2128 // get filename
2129 std::string filename = setup_filenames<map_t,engine_utils::toltecIO::toltec,engine_utils::toltecIO::psd>(dir_name);
2130
2131 // create file
2132 netCDF::NcFile fo(filename + ".nc", netCDF::NcFile::replace);
2133
2134 // loop through psd vector
2135 for (Eigen::Index i=0; i<mb->psds.size(); ++i) {
2136 // get name for extension layer
2137 std::string map_name = get_map_name(i);
2138
2139 // get the array for the given map
2140 Eigen::Index map_index = arrays_to_maps(i);
2141 // get the stokes parameter for the given map
2142 Eigen::Index stokes_index = maps_to_stokes(i);
2143
2144 auto array = calib.arrays[map_index];
2145 std::string name = toltec_io.array_name_map[array] + "_" + map_name + rtcproc.polarization.stokes_params[stokes_index];
2146
2147 // add dimensions
2148 netCDF::NcDim psd_dim = fo.addDim(name + "_nfreq",mb->psds[i].size());
2149 netCDF::NcDim pds_2d_row_dim = fo.addDim(name + "_rows",mb->psd_2ds[i].rows());
2150 netCDF::NcDim pds_2d_col_dim = fo.addDim(name + "_cols",mb->psd_2ds[i].cols());
2151
2152 std::vector<netCDF::NcDim> dims;
2153 dims.push_back(pds_2d_row_dim);
2154 dims.push_back(pds_2d_col_dim);
2155
2156 // psd
2157 netCDF::NcVar psd_v = fo.addVar(name + "_psd",netCDF::ncDouble, psd_dim);
2158 psd_v.putVar(mb->psds[i].data());
2159
2160 // psd freq
2161 netCDF::NcVar psd_freq_v = fo.addVar(name + "_psd_freq",netCDF::ncDouble, psd_dim);
2162 psd_freq_v.putVar(mb->psd_freqs[i].data());
2163
2164 // transpose 2d psd and freq
2165 Eigen::MatrixXd psd_2d_transposed = mb->psd_2ds[i].transpose();
2166 Eigen::MatrixXd psd_2d_freq_transposed = mb->psd_2d_freqs[i].transpose();
2167
2168 // 2d psd
2169 netCDF::NcVar psd_2d_v = fo.addVar(name + "_psd_2d",netCDF::ncDouble, dims);
2170 psd_2d_v.putVar(psd_2d_transposed.data());
2171
2172 // 2d psd freq
2173 netCDF::NcVar psd_2d_freq_v = fo.addVar(name + "_psd_2d_freq",netCDF::ncDouble, dims);
2174 psd_2d_freq_v.putVar(psd_2d_freq_transposed.data());
2175
2176 if (!mb->noise.empty()) {
2177 // add dimensions
2178 netCDF::NcDim noise_psd_dim = fo.addDim(name + "_noise_nfreq",mb->noise_psds[i].size());
2179 netCDF::NcDim noise_pds_2d_row_dim = fo.addDim(name + "_noise_rows",mb->noise_psd_2ds[i].rows());
2180 netCDF::NcDim noise_pds_2d_col_dim = fo.addDim(name + "_noise_cols",mb->noise_psd_2ds[i].cols());
2181
2182 std::vector<netCDF::NcDim> noise_dims;
2183 noise_dims.push_back(noise_pds_2d_row_dim);
2184 noise_dims.push_back(noise_pds_2d_col_dim);
2185
2186 // noise psd
2187 netCDF::NcVar noise_psd_v = fo.addVar(name + "_noise_psd",netCDF::ncDouble, noise_psd_dim);
2188 noise_psd_v.putVar(mb->noise_psds[i].data());
2189
2190 // noise psd freq
2191 netCDF::NcVar noise_psd_freq_v = fo.addVar(name + "_noise_psd_freq",netCDF::ncDouble, noise_psd_dim);
2192 noise_psd_freq_v.putVar(mb->noise_psd_freqs[i].data());
2193
2194 // transpose 2d noise psd and freq
2195 Eigen::MatrixXd noise_psd_2d_transposed = mb->noise_psd_2ds[i].transpose();
2196 Eigen::MatrixXd noise_psd_2d_freq_transposed = mb->noise_psd_2d_freqs[i].transpose();
2197
2198 // 2d noise psd
2199 netCDF::NcVar noise_psd_2d_v = fo.addVar(name + "_noise_psd_2d",netCDF::ncDouble, noise_dims);
2200 noise_psd_2d_v.putVar(noise_psd_2d_transposed.data());
2201
2202 // 2d noise psd freq
2203 netCDF::NcVar noise_psd_2d_freq_v = fo.addVar(name + "_noise_psd_2d_freq",netCDF::ncDouble, noise_dims);
2204 noise_psd_2d_freq_v.putVar(noise_psd_2d_freq_transposed.data());
2205 }
2206 }
2207 // close file
2208 fo.close();
2209}
2210
2211template <mapmaking::MapType map_t, class map_buffer_t>
2212void Engine::write_hist(map_buffer_t &mb, std::string dir_name) {
2213 std::string filename = setup_filenames<map_t,engine_utils::toltecIO::toltec,engine_utils::toltecIO::hist>(dir_name);
2214
2215 netCDF::NcFile fo(filename + ".nc", netCDF::NcFile::replace);
2216 netCDF::NcDim hist_bins_dim = fo.addDim("n_bins", mb->hist_n_bins);
2217
2218 // loop through stored histograms
2219 for (Eigen::Index i=0; i<mb->hists.size(); ++i) {
2220 // string to hold name
2221 // get name for extension layer
2222 std::string map_name = get_map_name(i);
2223
2224 // get the array for the given map
2225 Eigen::Index map_index = arrays_to_maps(i);
2226 // get the stokes parameter for the given map
2227 Eigen::Index stokes_index = maps_to_stokes(i);
2228
2229 // array index
2230 auto array = calib.arrays[map_index];
2231 std::string name = toltec_io.array_name_map[array] + "_" + map_name + rtcproc.polarization.stokes_params[stokes_index];
2232
2233 // histogram bins
2234 netCDF::NcVar hist_bins_v = fo.addVar(name + "_bins",netCDF::ncDouble, hist_bins_dim);
2235 hist_bins_v.putVar(mb->hist_bins[i].data());
2236
2237 // histogram
2238 netCDF::NcVar hist_v = fo.addVar(name + "_hist",netCDF::ncDouble, hist_bins_dim);
2239 hist_v.putVar(mb->hists[i].data());
2240
2241 if (!mb->noise.empty()) {
2242 // average noise histogram
2243 netCDF::NcVar hist_v = fo.addVar(name + "_noise_hist",netCDF::ncDouble, hist_bins_dim);
2244 hist_v.putVar(mb->noise_hists[i].data());
2245 }
2246 }
2247 // close file
2248 fo.close();
2249}
2250
2252 std::string path = obsnum_dir_name + "raw/";
2253 // if using tod subdir, put stats file in it
2254 if (tod_output_subdir_name!="null") {
2255 if (!fs::exists(fs::status(path + tod_output_subdir_name))) {
2256 fs::create_directories(path + tod_output_subdir_name);
2257 path = path + tod_output_subdir_name + "/";
2258 }
2259 }
2260 // create stats filename
2263 (path, redu_type, "", obsnum, telescope.sim_obs);
2264
2265 // det stats header
2266 std::map<std::string, std::string> det_stats_header_units {
2267 {"rms", omb.sig_unit},
2268 {"stddev",omb.sig_unit},
2269 {"median",omb.sig_unit},
2270 {"flagged_frac","N/A"},
2271 {"weights","1/(" + omb.sig_unit + ")^2"},
2272 };
2273 // group stats header
2274 std::map<std::string, std::string> grp_stats_header_units {
2275 {"median_weights", "1/(" + omb.sig_unit + ")^2"},
2276 };
2277
2278 netCDF::NcFile fo(stats_filename + ".nc", netCDF::NcFile::replace);
2279
2280 // add obsnum
2281 netCDF::NcVar obsnum_v = fo.addVar("obsnum",netCDF::ncInt);
2282 obsnum_v.putAtt("units","N/A");
2283 int obsnum_int = std::stoi(obsnum);
2284 obsnum_v.putVar(&obsnum_int);
2285
2286 // add dimensions
2287 netCDF::NcDim n_dets_dim = fo.addDim("n_dets", calib.n_dets);
2288 netCDF::NcDim n_arrays_dim = fo.addDim("n_arrays", calib.n_arrays);
2289 netCDF::NcDim n_chunks_dim = fo.addDim("n_chunks", telescope.scan_indices.cols());
2290
2291 std::vector<netCDF::NcDim> dims = {n_chunks_dim, n_dets_dim};
2292 std::vector<netCDF::NcDim> grp_dims = {n_chunks_dim, n_arrays_dim};
2293
2294 // add det stats
2295 for (const auto &stat: diagnostics.det_stats_header) {
2296 netCDF::NcVar stat_v = fo.addVar(stat,netCDF::ncDouble, dims);
2297 stat_v.putVar(diagnostics.stats[stat].data());
2298 stat_v.putAtt("units",det_stats_header_units[stat]);
2299 }
2300 // add group stats
2301 for (const auto &stat: diagnostics.grp_stats_header) {
2302 netCDF::NcVar stat_v = fo.addVar(stat,netCDF::ncDouble, grp_dims);
2303 stat_v.putVar(diagnostics.stats[stat].data());
2304 stat_v.putAtt("units",grp_stats_header_units[stat]);
2305 }
2306
2307 // add apt table
2308 for (auto const& x: calib.apt) {
2309 netCDF::NcVar apt_v = fo.addVar("apt_" + x.first,netCDF::ncDouble, n_dets_dim);
2310 apt_v.putVar(x.second.data());
2311 apt_v.putAtt("units",calib.apt_header_units[x.first]);
2312 }
2313
2314 // add adc
2315 if (!diagnostics.adc_snap_data.empty()) {
2316 netCDF::NcDim adc_snap_dim = fo.addDim("adcSnapDim", diagnostics.adc_snap_data[0].cols());
2317 netCDF::NcDim adc_snap_data_dim = fo.addDim("adcSnapDataDim", diagnostics.adc_snap_data[0].rows());
2318 std::vector<netCDF::NcDim> dims = {adc_snap_dim, adc_snap_data_dim};
2319 Eigen::Index i = 0;
2320 for (auto const& x: diagnostics.adc_snap_data) {
2321 netCDF::NcVar adc_snap_v = fo.addVar("toltec" + std::to_string(calib.nws(i)) + "_adc_snap_data",netCDF::ncDouble, dims);
2322 adc_snap_v.putVar(x.data());
2323 i++;
2324 }
2325 }
2326
2327 // add eigenvalues
2328 if (!diagnostics.evals.empty()) {
2329 netCDF::NcDim n_eigs_dim = fo.addDim("n_eigs",ptcproc.cleaner.n_calc);
2330 netCDF::NcDim n_eig_grp_dim = fo.addDim("n_eig_grp",diagnostics.evals[0][0].size());
2331
2332 std::vector<netCDF::NcDim> eval_dims = {n_eig_grp_dim, n_eigs_dim};
2333
2334 // loop through chunks
2335 for (const auto &[key, val]: diagnostics.evals) {
2336 // loop through cleaner gropuing
2337 for (Eigen::Index i=0; i<val.size(); ++i) {
2338
2339 netCDF::NcVar eval_v = fo.addVar("evals_" + ptcproc.cleaner.grouping[i] + "_" + std::to_string(i) +
2340 "_chunk_" + std::to_string(key), netCDF::ncDouble,eval_dims);
2341 std::vector<std::size_t> start_eig_index = {0, 0};
2342 std::vector<std::size_t> size = {1, TULA_SIZET(ptcproc.cleaner.n_calc)};
2343
2344 // loop through eigenvalues in current group
2345 for (const auto &evals: val[i]) {
2346 eval_v.putVar(start_eig_index,size,evals.data());
2347 start_eig_index[0] += 1;
2348 }
2349 }
2350 }
2351 }
2352 fo.close();
2353}
2354
2355template <mapmaking::MapType map_t, class map_buffer_t>
2356void Engine::run_wiener_filter(map_buffer_t &mb) {
2357 // pointer to map buffer
2358 mapmaking::MapBuffer* pmb = &mb;
2359 // pointer to data file fits vector
2360 std::vector<fitsIO<file_type_enum::write_fits, CCfits::ExtHDU*>>* f_io = nullptr;
2361 // pointer to noise file fits vector
2362 std::vector<fitsIO<file_type_enum::write_fits, CCfits::ExtHDU*>>* n_io = nullptr;
2363 // directory name
2364 std::string dir_name;
2365
2366 // filtered obs maps
2367 if constexpr (map_t == mapmaking::FilteredObs) {
2368 f_io = &filtered_fits_io_vec;
2370 dir_name = obsnum_dir_name + "filtered/";
2371 }
2372
2373 // filtered coadded maps
2374 else if constexpr (map_t == mapmaking::FilteredCoadd) {
2377 dir_name = coadd_dir_name + "filtered/";
2378 }
2379
2380 for (Eigen::Index i=0; i<f_io->size(); ++i) {
2381 // get the array for the given map
2382 // add primary hdu
2383 add_phdu(f_io, pmb, i);
2384
2385 // add primary hdu to noise maps
2386 if (!pmb->noise.empty()) {
2387 add_phdu(n_io, pmb, i);
2388 }
2389 }
2390
2391 // loop through maps and run wiener filter
2392 for (Eigen::Index i=0; i<n_maps; ++i) {
2393 // current array
2394 auto array = maps_to_arrays(i);
2395 // get file index
2396 auto map_index = arrays_to_maps(i);
2397 // init fwhm in pixels
2399 // make wiener filter template
2401 // run the filter for the current map
2403
2404 // filter noise maps
2405 if (run_noise) {
2406 tula::logging::progressbar pb(
2407 [&](const auto &msg) { logger->info("{}", msg); }, 100,
2408 "filtering noise");
2409
2410 for (Eigen::Index j=0; j<mb.n_noise; ++j) {
2411 wiener_filter.filter_noise(mb, i, j);
2412 pb.count(mb.n_noise, mb.n_noise / 100);
2413 }
2414
2416 logger->info("renormalizing errors");
2417 // get median error from weight maps
2418 mb.calc_median_err();
2419 // get median map rms from noise maps
2420 mb.calc_median_rms();
2421
2422 // get rescaled normalization factor
2423 auto noise_factor = (1./pow(mb.median_rms.array(),2.))*mb.median_err.array();
2424 // re-normalize weight map
2425 mb.weight[i].noalias() = mb.weight[i]*noise_factor(i);
2426
2427 logger->info("median rms {} ({})", static_cast<float>(mb.median_rms(i)), mb.sig_unit);
2428 }
2429 }
2430
2432 // only write if saving all iterations or on last iteration
2433 // write maps immediately after filtering due to computation time
2434 write_maps(f_io,n_io,pmb,i);
2435
2436 logger->info("file has been written to:");
2437 logger->info("{}.fits",f_io->at(map_index).filepath);
2438
2439 // explicitly destroy the fits file after we're done with it
2440 bool close_file = true;
2443 close_file = false;
2444 }
2445 }
2446 // check if we're moving onto a new file
2447 if (i<n_maps-1) {
2448 if (arrays_to_maps(i+1) > arrays_to_maps(i) && close_file) {
2449 f_io->at(map_index).pfits->destroy();
2450 }
2451 }
2452 }
2453 }
2454
2456 // clear fits file vectors to ensure its closed.
2457 f_io->clear();
2458 n_io->clear();
2459 }
2460}
2461
2462template <mapmaking::MapType map_t, class map_buffer_t>
2463void Engine::find_sources(map_buffer_t &mb) {
2464 // clear all source vectors
2465 mb.n_sources.clear();
2466 mb.row_source_locs.clear();
2467 mb.col_source_locs.clear();
2468 // loop through maps
2469 for (Eigen::Index i=0; i<n_maps; ++i) {
2470 // update source vectors
2471 mb.n_sources.push_back(0);
2472 mb.row_source_locs.push_back(Eigen::VectorXi::Ones(1));
2473 mb.col_source_locs.push_back(Eigen::VectorXi::Ones(1));
2474
2475 // default value of -99 to keep size of vectors same as map vector
2476 mb.row_source_locs.back()*=-99;
2477 mb.col_source_locs.back()*=-99;
2478
2479 // run source finder
2480 auto sources_found = mb.find_sources(i);
2481
2482 // number of sources found for current map
2483 if (sources_found) {
2484 logger->info("{} source(s) found", mb.n_sources.back());
2485 }
2486 else {
2487 logger->info("no sources found");
2488 }
2489 }
2490
2491 // count up the total number of sources
2492 Eigen::Index n_sources = 0;
2493 for (const auto &sources: mb.n_sources) {
2494 n_sources += sources;
2495 }
2496
2497 // matrix to store source parameters
2498 mb.source_params.setZero(n_sources,map_fitter.n_params);
2499 mb.source_perror.setZero(n_sources,map_fitter.n_params);
2500
2501 // keep track of row in total source count
2502 Eigen::Index k = 0;
2503
2504 // now loop through and fit the sources
2505 for (Eigen::Index i=0; i<n_maps; ++i) {
2506 // skip map if no sources found
2507 if (mb.n_sources[i] > 0) {
2508 // current array
2509 auto array = maps_to_arrays(i);
2510 // init fwhm in pixels
2511 auto init_fwhm = toltec_io.array_fwhm_arcsec[array]*ASEC_TO_RAD/mb.pixel_size_rad;
2512
2513 // placeholder vectors for grppi map
2514 std::vector<int> source_in_vec, source_out_vec;
2515
2516 source_in_vec.resize(mb.n_sources[i]);
2517 std::iota(source_in_vec.begin(), source_in_vec.end(), 0);
2518 source_out_vec.resize(mb.n_sources[i]);
2519
2520 // loop through sources and fit them
2521 grppi::map(tula::grppi_utils::dyn_ex(parallel_policy), source_in_vec, source_out_vec, [&](auto j) {
2522 // update source rows and cols
2523 double init_row = mb.row_source_locs[i](j);
2524 double init_col = mb.col_source_locs[i](j);
2525
2526 // fit source
2527 auto [params, perrors, good_fit] =
2529 init_fwhm, init_row, init_col);
2530 if (good_fit) {
2531 // rescale fit params from pixel to on-sky units
2532 params(1) = RAD_TO_ASEC*mb.pixel_size_rad*(params(1) - (mb.n_cols)/2);
2533 params(2) = RAD_TO_ASEC*mb.pixel_size_rad*(params(2) - (mb.n_rows)/2);
2534 params(3) = RAD_TO_ASEC*STD_TO_FWHM*mb.pixel_size_rad*(params(3));
2535 params(4) = RAD_TO_ASEC*STD_TO_FWHM*mb.pixel_size_rad*(params(4));
2536
2537 // rescale fit errors from pixel to on-sky units
2538 perrors(1) = RAD_TO_ASEC*mb.pixel_size_rad*(perrors(1));
2539 perrors(2) = RAD_TO_ASEC*mb.pixel_size_rad*(perrors(2));
2540 perrors(3) = RAD_TO_ASEC*STD_TO_FWHM*mb.pixel_size_rad*(perrors(3));
2541 perrors(4) = RAD_TO_ASEC*STD_TO_FWHM*mb.pixel_size_rad*(perrors(4));
2542
2543 // if in radec calculate absolute pointing
2544 if (telescope.pixel_axes=="radec") {
2545 Eigen::VectorXd lat(1), lon(1);
2546 lat << params(2)*ASEC_TO_RAD;
2547 lon << params(1)*ASEC_TO_RAD;
2548
2549 auto [adec, ara] = engine_utils::tangent_to_abs(lat, lon, mb.wcs.crval[0]*DEG_TO_RAD, mb.wcs.crval[1]*DEG_TO_RAD);
2550
2551 params(1) = ara(0)*RAD_TO_DEG;
2552 params(2) = adec(0)*RAD_TO_DEG;
2553
2554 perrors(1) = perrors(1)*ASEC_TO_DEG;
2555 perrors(2) = perrors(2)*ASEC_TO_DEG;
2556 }
2557
2558 // add source params and errors to table
2559 mb.source_params.row(k+j) = params;
2560 mb.source_perror.row(k+j) = perrors;
2561 }
2562 return 0;
2563 });
2564
2565 // update row
2566 k += mb.n_sources[i];
2567 }
2568 }
2569}
2570
2571template <mapmaking::MapType map_t, class map_buffer_t>
2572void Engine::write_sources(map_buffer_t &mb, std::string dir_name) {
2573 // get filenmame for source table
2574 std::string source_filename = setup_filenames<map_t,engine_utils::toltecIO::source,
2575 engine_utils::toltecIO::map>(dir_name);
2576
2577 // source header information
2578 std::vector<std::string> source_header = {
2579 "array",
2580 "amp",
2581 "amp_err",
2582 "x_t",
2583 "x_t_err",
2584 "y_t",
2585 "y_t_err",
2586 "a_fwhm",
2587 "a_fwhm_err",
2588 "b_fwhm",
2589 "b_fwhm_err",
2590 "angle",
2591 "angle_err",
2592 "sig2noise"
2593 };
2594
2595 // units for fitted parameter centroids
2596 std::string pos_units = (telescope.pixel_axes == "radec") ? "deg" : "arcsec";
2597
2598 // units for source header
2599 std::map<std::string,std::string> source_header_units = {
2600 {"array","N/A"},
2601 {"amp", mb->sig_unit},
2602 {"amp_err", mb->sig_unit},
2603 {"x_t", pos_units},
2604 {"x_t_err", pos_units},
2605 {"y_t", pos_units},
2606 {"y_t_err", pos_units},
2607 {"a_fwhm", "arcsec"},
2608 {"a_fwhm_err", "arcsec"},
2609 {"b_fwhm", "arcsec"},
2610 {"b_fwhm_err", "arcsec"},
2611 {"angle", "rad"},
2612 {"angle_err", "rad"},
2613 {"sig2noise", "N/A"}
2614 };
2615
2616 // meta information for source table
2617 YAML::Node source_meta;
2618
2619 // add obsnums
2620 for (Eigen::Index i=0; i<mb->obsnums.size(); ++i) {
2621 // add obsnum to meta data
2622 source_meta["obsnum" + std::to_string(i)] = mb->obsnums[i];
2623 }
2624
2625 // add source name
2626 source_meta["Source"] = telescope.source_name;
2627
2628 // add date of file creation
2629 source_meta["creation_date"] = engine_utils::current_date_time();
2630
2631 // add observation date
2632 source_meta["date"] = date_obs.back();
2633
2634
2635 // populate source meta information
2636 for (const auto &[key,val]: source_header_units) {
2637 source_meta[key].push_back("units: " + val);
2638 // description from apt
2639 auto description = calib.apt_header_description[key];
2640 source_meta[key].push_back(description);
2641 }
2642
2643 // count up the total number of sources
2644 Eigen::Index n_sources = 0;
2645 for (const auto &sources: mb->n_sources) {
2646 n_sources += sources;
2647 }
2648
2649 // matrix to hold source information (floats for readability)
2650 Eigen::MatrixXf source_table(n_sources, 2*map_fitter.n_params + 2);
2651
2652 // loop through params and add arrays
2653 Eigen::Index k=0;
2654 for (Eigen::Index i=0; i<mb->n_sources.size(); ++i) {
2655 if (mb->n_sources[i]!=0) {
2656 // calculate map standard deviation
2657 double map_std_dev = engine_utils::calc_std_dev(mb->signal[i]);
2658
2659 for (Eigen::Index j=0; j<mb->n_sources[i]; ++j) {
2660 source_table(k,0) = maps_to_arrays(i);
2661 // set signal to noise
2662 source_table(k,2*map_fitter.n_params + 1) = mb->source_params(k,0)/map_std_dev;
2663
2664 k++;
2665 }
2666 }
2667 }
2668
2669 // populate source table
2670 Eigen::Index j = 0;
2671 for (Eigen::Index i=1; i<2*map_fitter.n_params; i=i+2) {
2672 source_table.col(i) = mb->source_params.col(j).template cast <float> ();
2673 source_table.col(i+1) = mb->source_perror.col(j).template cast <float> ();
2674 j++;
2675 }
2676
2677 // write source table
2678 to_ecsv_from_matrix(source_filename, source_table, source_header, source_meta);
2679}
Definition engine.h:155
Eigen::Index hwpr_start_indices
Definition engine.h:208
void obsnum_setup()
Definition engine.h:353
std::string obsnum_dir_name
Definition engine.h:184
void write_stats()
Definition engine.h:2251
std::string parallel_policy
Definition engine.h:196
void write_map_summary(map_buffer_t &)
Definition engine.h:1584
std::string redu_dir_name
Definition engine.h:178
std::map< std::string, Eigen::VectorXd > pointing_offsets_arcsec
Definition engine.h:241
std::vector< fitsIO< file_type_enum::write_fits, CCfits::ExtHDU * > > filtered_coadd_fits_io_vec
Definition engine.h:251
std::vector< Eigen::Index > start_indices
Definition engine.h:205
void create_tod_files()
Definition engine.h:1276
std::vector< Eigen::VectorXd > nw_times
Definition engine.h:166
void write_psd(map_buffer_t &, std::string)
Definition engine.h:2127
int redu_dir_num
Definition engine.h:181
auto get_map_name(int)
Definition engine.h:1702
std::vector< fitsIO< file_type_enum::write_fits, CCfits::ExtHDU * > > noise_fits_io_vec
Definition engine.h:246
std::vector< fitsIO< file_type_enum::write_fits, CCfits::ExtHDU * > > coadd_noise_fits_io_vec
Definition engine.h:250
void get_beammap_config(CT &)
Definition engine.h:667
key_vec_t invalid_keys
Definition engine.h:190
std::string tod_output_type
Definition engine.h:223
auto setup_filenames(std::string dir_name)
Definition engine.h:1674
std::string tod_output_subdir_name
Definition engine.h:223
void create_obs_map_files()
Definition engine.h:991
key_vec_t missing_keys
Definition engine.h:190
std::string coadd_dir_name
Definition engine.h:184
std::vector< fitsIO< file_type_enum::write_fits, CCfits::ExtHDU * > > filtered_noise_fits_io_vec
Definition engine.h:247
void get_mapmaking_config(CT &)
Definition engine.h:559
Eigen::VectorXI maps_to_arrays
Definition engine.h:232
int fruit_iter
Definition engine.h:238
bool verbose_mode
Definition engine.h:172
int n_maps
Definition engine.h:229
Eigen::VectorXI arrays_to_maps
Definition engine.h:232
std::vector< fitsIO< file_type_enum::write_fits, CCfits::ExtHDU * > > filtered_fits_io_vec
Definition engine.h:247
void run_wiener_filter(map_buffer_t &)
Definition engine.h:2356
std::string output_dir
Definition engine.h:178
int n_scans_done
Definition engine.h:199
void add_tod_header()
Definition engine.h:1054
bool write_filtered_maps_partial
Definition engine.h:220
Eigen::Index hwpr_end_indices
Definition engine.h:208
std::vector< Eigen::Index > end_indices
Definition engine.h:205
Eigen::VectorXd t_common
Definition engine.h:164
void write_hist(map_buffer_t &, std::string)
Definition engine.h:2212
std::vector< fitsIO< file_type_enum::write_fits, CCfits::ExtHDU * > > filtered_coadd_noise_fits_io_vec
Definition engine.h:251
void add_phdu(fits_io_type &, map_buffer_t &, Eigen::Index)
Definition engine.h:1737
std::string redu_type
Definition engine.h:214
std::map< std::string, int > gaps
Definition engine.h:175
std::vector< Eigen::VectorXi > masks
Definition engine.h:165
std::string map_grouping
Definition engine.h:226
std::shared_ptr< spdlog::logger > logger
Definition engine.h:161
void get_citlali_config(CT &)
Definition engine.h:772
std::map< std::string, std::string > tod_filename
Definition engine.h:187
void write_sources(map_buffer_t &, std::string)
Definition engine.h:2572
void get_timestream_config(CT &)
Definition engine.h:497
std::string tod_type
Definition engine.h:211
void cli_summary()
Definition engine.h:1463
Eigen::VectorXI maps_to_stokes
Definition engine.h:235
void get_astrometry_config(CT &)
Definition engine.h:965
void get_photometry_config(CT &)
Definition engine.h:931
void get_rtc_config(CT &)
Definition engine.h:466
std::string obsnum
Definition engine.h:217
std::vector< fitsIO< file_type_enum::write_fits, CCfits::ExtHDU * > > coadd_fits_io_vec
Definition engine.h:250
void find_sources(map_buffer_t &)
Definition engine.h:2463
std::map< std::string, double > interface_sync_offset
Definition engine.h:202
std::vector< std::string > date_obs
Definition engine.h:169
void write_chunk_summary(TCData< tc_t, Eigen::MatrixXd > &)
Definition engine.h:1529
Eigen::ArrayXd pointing_offsets_modified_julian_date
Definition engine.h:243
void get_ptc_config(CT &)
Definition engine.h:486
std::vector< fitsIO< file_type_enum::write_fits, CCfits::ExtHDU * > > fits_io_vec
Definition engine.h:246
int n_threads
Definition engine.h:193
std::string map_method
Definition engine.h:226
std::vector< std::vector< std::string > > key_vec_t
Definition engine.h:158
void get_map_filter_config(CT &)
Definition engine.h:737
void write_maps(fits_io_type &, fits_io_type &, map_buffer_t &, Eigen::Index)
Definition engine.h:2030
Definition calib.h:14
Eigen::VectorXI fg
Definition calib.h:29
bool run_hwpr
Definition calib.h:41
std::string ignore_hwpr
Definition calib.h:32
Eigen::Index n_arrays
Definition calib.h:53
std::map< Eigen::Index, double > array_pas
Definition calib.h:62
std::map< Eigen::Index, double > array_beam_areas
Definition calib.h:65
Eigen::VectorXI nws
Definition calib.h:44
Eigen::VectorXI arrays
Definition calib.h:44
std::map< std::string, std::string > apt_header_units
Definition calib.h:116
std::string apt_filepath
Definition calib.h:20
std::map< std::string, Eigen::VectorXd > apt
Definition calib.h:22
std::map< std::string, std::string > apt_header_description
Definition calib.h:151
std::map< Eigen::Index, std::tuple< double, double > > array_fwhms
Definition calib.h:59
Eigen::Index n_dets
Definition calib.h:47
Definition diagnostics.h:14
std::map< Eigen::Index, std::vector< std::vector< Eigen::VectorXd > > > evals
Definition diagnostics.h:45
bool write_evals
Definition diagnostics.h:25
std::vector< std::string > grp_stats_header
Definition diagnostics.h:40
std::map< std::string, Eigen::MatrixXd > stats
Definition diagnostics.h:16
std::vector< std::string > det_stats_header
Definition diagnostics.h:31
std::vector< Eigen::Matrix< short, Eigen::Dynamic, Eigen::Dynamic > > adc_snap_data
Definition diagnostics.h:19
Definition telescope.h:12
bool sim_obs
Definition telescope.h:20
std::map< std::string, Eigen::VectorXd > tel_data
Definition telescope.h:54
double fsmp
Definition telescope.h:42
std::string chunk_mode
Definition telescope.h:30
std::string obs_goal
Definition telescope.h:23
std::string source_name
Definition telescope.h:23
double tau_225_GHz
Definition telescope.h:45
Eigen::MatrixXI scan_indices
Definition telescope.h:51
std::string project_id
Definition telescope.h:23
bool force_chunk
Definition telescope.h:36
std::string pixel_axes
Definition telescope.h:59
std::map< std::string, Eigen::VectorXd > tel_header
Definition telescope.h:56
double chunking_value
Definition telescope.h:33
Eigen::Index inner_scans_chunk
Definition telescope.h:48
Definition fitting.h:13
double bounding_box_pix
Definition fitting.h:26
auto fit_to_gaussian(Eigen::DenseBase< Derived > &, Eigen::DenseBase< Derived > &, double, double, double)
Definition fitting.h:172
int n_params
Definition fitting.h:23
double fwhm_low
Definition fitting.h:38
double flux_low
Definition fitting.h:33
double fitting_region_pix
Definition fitting.h:27
bool fit_angle
Definition fitting.h:43
@ pointing
Definition fitting.h:18
double fwhm_high
Definition fitting.h:40
double flux_high
Definition fitting.h:35
Eigen::VectorXd fwhm_limits
Definition fitting.h:30
Eigen::VectorXd flux_limits
Definition fitting.h:30
Definition toltec_io.h:10
ProdType
Definition toltec_io.h:21
@ map
Definition toltec_io.h:22
@ stats
Definition toltec_io.h:28
@ ptc_timestream
Definition toltec_io.h:27
@ rtc_timestream
Definition toltec_io.h:26
@ noise
Definition toltec_io.h:23
std::map< Eigen::Index, std::string > array_name_map
Definition toltec_io.h:54
std::map< Eigen::Index, double > array_freq_map
Definition toltec_io.h:68
std::map< Eigen::Index, double > array_wavelength_map
Definition toltec_io.h:61
DataType
Definition toltec_io.h:13
@ source
Definition toltec_io.h:17
@ toltec
Definition toltec_io.h:14
@ filtered
Definition toltec_io.h:34
@ raw
Definition toltec_io.h:33
std::string create_filename(const std::string, const std::string, std::string, std::string, const bool)
Definition toltec_io.h:86
std::map< Eigen::Index, double > array_fwhm_arcsec
Definition toltec_io.h:75
Definition fits_io.h:10
Definition jinc_mm.h:26
void calculate_jinc_splines()
Definition jinc_mm.h:152
std::string mode
Definition jinc_mm.h:56
std::map< Eigen::Index, Eigen::VectorXd > shape_params
Definition jinc_mm.h:68
double r_max
Definition jinc_mm.h:62
void allocate_jinc_matrix(double)
Definition jinc_mm.h:118
std::string parallel_policy
Definition jinc_mm.h:53
bool run_polarization
Definition jinc_mm.h:50
Definition ml_mm.h:22
double tolerance
Definition ml_mm.h:27
int max_iterations
Definition ml_mm.h:28
Definition map.h:45
WCS wcs
Definition map.h:48
double exposure_time
Definition map.h:75
std::string sig_unit
Definition map.h:73
double source_window_rad
Definition map.h:128
Eigen::Index n_rows
Definition map.h:65
std::string parallel_policy
Definition map.h:59
Eigen::Index n_cols
Definition map.h:65
std::vector< Eigen::MatrixXd > coverage
Definition map.h:78
std::vector< Eigen::MatrixXd > signal
Definition map.h:78
std::string source_finder_mode
Definition map.h:123
std::vector< Eigen::MatrixXd > weight
Definition map.h:78
bool randomize_dets
Definition map.h:87
std::vector< Eigen::MatrixXd > kernel
Definition map.h:78
Eigen::Index n_noise
Definition map.h:67
void get_config(tula::config::YamlConfig &, std::vector< std::vector< std::string > > &, std::vector< std::vector< std::string > > &, std::string, std::string)
Definition map.cpp:15
double source_sigma
Definition map.h:126
std::vector< Eigen::Tensor< double, 3 > > noise
Definition map.h:81
double pixel_size_rad
Definition map.h:69
Definition naive_mm.h:20
bool run_polarization
Definition naive_mm.h:52
Definition wiener_filter.h:23
void filter_noise(MB &mb, const int, const int)
Definition wiener_filter.h:939
std::string template_type
Definition wiener_filter.h:29
std::string parallel_policy
Definition wiener_filter.h:59
std::map< std::string, double > template_fwhm_rad
Definition wiener_filter.h:50
std::string filter_type
Definition wiener_filter.h:29
bool normalize_error
Definition wiener_filter.h:31
void filter_maps(MB &, const int)
Definition wiener_filter.h:870
double init_fwhm
Definition wiener_filter.h:47
engine_utils::mapFitter map_fitter
Definition wiener_filter.h:69
void get_config(config_t &, std::vector< std::vector< std::string > > &, std::vector< std::vector< std::string > > &)
Definition wiener_filter.h:126
void make_template(MB &, CD &c, const double, const int)
Definition wiener_filter.h:694
bool run_lowpass
Definition wiener_filter.h:35
auto calc_tau(Eigen::DenseBase< Derived > &, double)
Definition calibrate.h:136
void setup(double tau_225_zenith)
Definition calibrate.h:31
std::string extinction_model
Definition calibrate.h:18
std::map< Eigen::Index, Eigen::VectorXI > n_eig_to_cut
Definition clean.h:27
int n_calc
Definition clean.h:30
std::vector< std::string > grouping
Definition clean.h:36
double fsmp
Definition despike.h:22
void make_filter(double)
Definition filter.h:34
Eigen::Index n_terms
Definition filter.h:19
std::string type
Definition kernel.h:15
std::string map_grouping
Definition kernel.h:28
double fwhm_rad
Definition kernel.h:22
void setup(Eigen::Index)
Definition kernel.h:50
Definition ptcproc.h:19
void get_config(config_t &, std::vector< std::vector< std::string > > &, std::vector< std::vector< std::string > > &)
Definition ptcproc.h:62
bool run_clean
Definition ptcproc.h:22
timestream::Cleaner cleaner
Definition ptcproc.h:31
double upper_weight_factor
Definition ptcproc.h:26
double lower_weight_factor
Definition ptcproc.h:26
std::string weighting_type
Definition ptcproc.h:28
double med_weight_factor
Definition ptcproc.h:24
std::map< int, std::string > stokes_params
Definition polarization.h:15
Definition rtcproc.h:18
bool run_polarization
Definition rtcproc.h:23
bool run_kernel
Definition rtcproc.h:24
void get_config(config_t &, std::vector< std::vector< std::string > > &, std::vector< std::vector< std::string > > &)
Definition rtcproc.h:71
bool run_calibrate
Definition rtcproc.h:28
bool run_downsample
Definition rtcproc.h:27
timestream::Filter filter
Definition rtcproc.h:35
bool run_tod_filter
Definition rtcproc.h:26
bool run_despike
Definition rtcproc.h:25
bool run_extinction
Definition rtcproc.h:29
timestream::Calibration calibration
Definition rtcproc.h:37
timestream::Despiker despiker
Definition rtcproc.h:34
timestream::Kernel kernel
Definition rtcproc.h:33
timestream::Polarization polarization
Definition rtcproc.h:32
std::string fruit_loops_path
Definition timestream.h:281
double upper_inv_var_factor
Definition timestream.h:302
double fruit_loops_sig2noise
Definition timestream.h:289
bool run_tod_output
Definition timestream.h:276
bool write_evals
Definition timestream.h:276
int fruit_loops_iters
Definition timestream.h:287
bool run_fruit_loops
Definition timestream.h:279
Eigen::VectorXd fruit_loops_flux
Definition timestream.h:291
double lower_inv_var_factor
Definition timestream.h:302
std::string fruit_loops_type
Definition timestream.h:285
void get_config_value(config_t config, param_t &param, key_vec_t &missing_keys, key_vec_t &invalid_keys, option_t option, std::vector< param_t > allowed={}, std::vector< param_t > min_val={}, std::vector< param_t > max_val={})
Definition config.h:58
#define DEG_TO_RAD
Definition constants.h:27
#define MJY_SR_TO_mJY_ASEC
Definition constants.h:51
#define STD_TO_FWHM
Definition constants.h:45
#define FWHM_TO_STD
Definition constants.h:48
#define RAD_TO_ASEC
Definition constants.h:36
#define ASEC_TO_RAD
Definition constants.h:33
#define RAD_TO_DEG
Definition constants.h:30
#define ASEC_TO_DEG
Definition constants.h:24
constexpr auto pi
Definition constants.h:6
void to_ecsv_from_matrix(std::string filepath, Eigen::DenseBase< Derived > &table, std::vector< std::string > header, YAML::Node meta)
Definition ecsv_io.h:50
static const double mJy_beam_to_uK(const double flux_mjy_beam, const double freq_Hz, const double fwhm)
Definition utils.h:77
auto tangent_to_abs(Eigen::DenseBase< Derived > &lat, Eigen::DenseBase< Derived > &lon, const double cra, const double cdec)
Definition pointing.h:88
static const int get_phys_memory()
Definition utils.h:56
static const std::string current_date_time()
Definition utils.h:91
auto calc_std_dev(Eigen::DenseBase< DerivedA > &data)
Definition utils.h:336
MapType
Definition map.h:17
@ FilteredObs
Definition map.h:19
@ FilteredCoadd
Definition map.h:21
@ RawObs
Definition map.h:18
@ RawCoadd
Definition map.h:20
void add_netcdf_var(netCDF::NcFile &fo, std::string name, T data)
Definition netcdf_io.h:11
Definition engine.h:115
bool beammap_subtract_reference
Definition engine.h:133
std::map< std::string, double > lower_sig2noise
Definition engine.h:148
std::map< std::string, double > beammap_fluxes_mJy_beam
Definition engine.h:123
Eigen::Index beammap_reference_det
Definition engine.h:136
std::map< std::string, double > upper_fwhm_arcsec
Definition engine.h:148
std::map< std::string, double > max_dist_arcsec
Definition engine.h:149
double lower_sens_factor
Definition engine.h:152
std::map< std::string, double > lower_fwhm_arcsec
Definition engine.h:148
double beammap_dec_rad
Definition engine.h:120
Eigen::VectorXd sens_psd_limits_Hz
Definition engine.h:145
int beammap_tod_output_iter
Definition engine.h:142
std::string beammap_source_name
Definition engine.h:117
std::map< std::string, double > beammap_err_MJy_Sr
Definition engine.h:124
double beammap_ra_rad
Definition engine.h:120
std::map< std::string, double > beammap_fluxes_MJy_Sr
Definition engine.h:124
int beammap_iter_max
Definition engine.h:127
bool beammap_derotate
Definition engine.h:139
std::map< std::string, double > upper_sig2noise
Definition engine.h:149
double beammap_iter_tolerance
Definition engine.h:130
std::map< std::string, double > beammap_err_mJy_beam
Definition engine.h:123
double upper_sens_factor
Definition engine.h:152
std::vector< float > crval
Definition map.h:36
Definition engine.h:93
mapmaking::NaiveMapmaker naive_mm
Definition engine.h:109
engine::Telescope telescope
Definition engine.h:96
engine_utils::mapFitter map_fitter
Definition engine.h:99
mapmaking::MLMapmaker ml_mm
Definition engine.h:111
mapmaking::WienerFilter wiener_filter
Definition engine.h:112
engine::Diagnostics diagnostics
Definition engine.h:98
engine_utils::toltecIO toltec_io
Definition engine.h:97
timestream::PTCProc ptcproc
Definition engine.h:105
mapmaking::JincMapmaker jinc_mm
Definition engine.h:110
mapmaking::MapBuffer omb
Definition engine.h:108
mapmaking::MapBuffer cmb
Definition engine.h:108
timestream::RTCProc rtcproc
Definition engine.h:102
engine::Calib calib
Definition engine.h:95
Definition engine.h:71
bool run_noise
Definition engine.h:86
bool use_subdir
Definition engine.h:75
bool run_coadd
Definition engine.h:85
bool run_map_filter
Definition engine.h:87
bool interp_over_gaps
Definition engine.h:73
bool run_source_finder
Definition engine.h:90
bool run_tod_output
Definition engine.h:81
bool run_mapmaking
Definition engine.h:84
bool run_tod
Definition engine.h:78
TC data class.
Definition timestream.h:55