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