Citlali
Loading...
Searching...
No Matches
beammap.h
Go to the documentation of this file.
1#pragma once
2
3#include <map>
4#include <vector>
5#include <string>
6
8
12
13// selects the type of TCData
14using timestream::TCDataKind;
15
16class Beammap: public Engine {
17public:
18 // parallel policies for each section
20
21 // vector to store each scan's PTCData
22 std::vector<TCData<TCDataKind::PTC,Eigen::MatrixXd>> ptcs0, ptcs;
23
24 // copy of obs map buffer for map iteration
26
27 // vector to store each scan's calib class
28 std::vector<engine::Calib> calib_scans0, calib_scans;
29
30 // beammap iteration parameters
31 Eigen::Index current_iter;
32
33 // vector for convergence check
34 Eigen::Matrix<bool, Eigen::Dynamic, 1> converged;
35
36 // vector to record convergence iteration
37 Eigen::Vector<int, Eigen::Dynamic> converge_iter;
38
39 // previous iteration fit parameters
40 Eigen::MatrixXd p0, perror0;
41
42 // current iteration fit parameters
43 Eigen::MatrixXd params, perrors;
44
45 // reference detector
46 Eigen::Index beammap_reference_det_found = -99;
47
48 // bitwise flags
49 enum AptFlags {
50 Good = 0,
51 BadFit = 1 << 0,
52 AzFWHM = 1 << 1,
53 ElFWHM = 1 << 2,
54 Sig2Noise = 1 << 3,
55 Sens = 1 << 4,
56 Position = 1 << 5,
57 };
58
59 // holds bitwise flags
60 Eigen::Matrix<uint16_t,Eigen::Dynamic,1> flag2;
61
62 // good fits
63 Eigen::Matrix<bool, Eigen::Dynamic, 1> good_fits;
64
65 // placeholder vectors for grppi maps
66 std::vector<int> scan_in_vec, scan_out_vec;
67 std::vector<int> det_in_vec, det_out_vec;
68
69 // initial setup for each obs
70 void setup();
71
72 // timestream grppi pipeline
73 template <class KidsProc, class RawObs>
74 void timestream_pipeline(KidsProc &, RawObs &);
75
76 // run the raw time chunk processing
77 template <class KidsProc>
78 auto run_timestream(KidsProc &);
79
80 // run the loop pipeline
81 void loop_pipeline();
82
83 // run the iterative stage
84 void run_loop();
85
86 // flag detectors
87 void set_apt_flags();
88
89 // derotate apt and subtract reference detector
90 void process_apt();
91
92 // main pipeline process
93 template <class KidsProc, class RawObs>
94 void pipeline(KidsProc &, RawObs &);
95
96 // output files
97 template <mapmaking::MapType map_type>
98 void output();
99};
100
102 // assign parallel policies
104
105 // run obsnum setup
106 obsnum_setup();
107
108 // create kids tone apt row
109 calib.apt["kids_tone"].resize(calib.n_dets);
110
111 Eigen::Index j = 0;
112 // set kids tone (det number on network)
113 calib.apt["kids_tone"](0) = 0;
114 for (Eigen::Index i=1; i<calib.n_dets; ++i) {
115 if (calib.apt["nw"](i) > calib.apt["nw"](i-1)) {
116 j = 0;
117 }
118 else {
119 j++;
120 }
121
122 calib.apt["kids_tone"](i) = j;
123 }
124
125 // add kids tone to apt header
126 calib.apt_header_keys.push_back("kids_tone");
127 calib.apt_header_units["kids_tone"] = "N/A";
128
129 // resize the PTCData vector to number of scans
130 ptcs0.resize(telescope.scan_indices.cols());
131
132 // resize the calib vector to number of scans
133 calib_scans0.resize(telescope.scan_indices.cols());
134
135 // resize the initial fit matrix
136 p0.setZero(n_maps, map_fitter.n_params);
137 // resize the initial fit error matrix
139 // resize the current fit matrix
142
143 // resize good fits
144 good_fits.setZero(n_maps);
145
146 // initially all detectors are unconverged
147 converged.setZero(n_maps);
148 // convergence iteration
149 converge_iter.resize(n_maps);
150 converge_iter.setConstant(1);
151 // set the initial iteration
152 current_iter = 0;
153
154 /* update apt table meta data */
155 calib.apt_meta.reset();
156
157 // add obsnum to meta data
158 calib.apt_meta["obsnum"] = obsnum;
159
160 // add source name
162
163 // add project id to meta data
164 calib.apt_meta["project_id"] = telescope.project_id;
165
166 // add input source flux
167 for (const auto &beammap_flux: beammap_fluxes_mJy_beam) {
168 auto key = beammap_flux.first + "_flux";
169 calib.apt_meta[key].push_back(beammap_flux.second);
170 calib.apt_meta[key].push_back("units: mJy/beam");
171 calib.apt_meta[key].push_back(beammap_flux.first + " flux density");
172 }
173
174 // add date of file creation
176
177 // add observation date
178 calib.apt_meta["date"] = date_obs.back();
179
180 // mean Modified Julian Date
182
183 // reference frame
184 calib.apt_meta["Radesys"] = telescope.pixel_axes;
185
186 // add mean tau to apt meta
188 Eigen::VectorXd tau_el(1);
189 tau_el << telescope.tel_data["TelElAct"].mean();
190 auto tau_freq = rtcproc.calibration.calc_tau(tau_el, telescope.tau_225_GHz);
191
192 Eigen::Index i = 0;
193 for (auto const& [key, val] : tau_freq) {
194 calib.apt_meta[toltec_io.array_name_map[calib.arrays(i)]+"_tau"] = val[0];
195 i++;
196 }
197 }
198 else {
199 for (Eigen::Index i=0; i<calib.arrays.size(); ++i) {
201 }
202 }
203
204 // add apt header keys
205 for (const auto &[param,unit]: calib.apt_header_units) {
206 calib.apt_meta[param].push_back("units: " + unit);
207 }
208 // add apt header descriptions
209 for (const auto &[param,description]: calib.apt_header_description) {
210 calib.apt_meta[param].push_back(description);
211 }
212
213 // kids tone
214 calib.apt_meta["kids_tone"].push_back("units: N/A");
215 calib.apt_meta["kids_tone"].push_back("index of tone in network");
216
217 // bitwise flag
218 calib.apt_meta["flag2"].push_back("units: N/A");
219 calib.apt_meta["flag2"].push_back("bitwise flag");
220 calib.apt_meta["flag2"].push_back("Good=0");
221 calib.apt_meta["flag2"].push_back("BadFit=1");
222 calib.apt_meta["flag2"].push_back("AzFWHM=2");
223 calib.apt_meta["flag2"].push_back("ElFWHM=3");
224 calib.apt_meta["flag2"].push_back("Sig2Noise=4");
225 calib.apt_meta["flag2"].push_back("Sens=5");
226 calib.apt_meta["flag2"].push_back("Position=6");
227
228 // add array mapping
229 for (const auto &[arr_index,arr_name]: toltec_io.array_name_map) {
230 calib.apt_meta["array_order"].push_back(std::to_string(arr_index) + ": " + arr_name);
231 }
232
233 calib.apt_header_units["flag2"] = "N/A";
234 calib.apt_header_keys.push_back("flag2");
235
236 // is the detector rotated?
237 calib.apt_meta["is_derotated"] = beammap_derotate;
238 // was a reference detector subtracted?
239 calib.apt_meta["reference_detector_subtracted"] = beammap_subtract_reference;
240 // reference detector
241 calib.apt_meta["reference_det"] = beammap_reference_det_found;
242}
243
244template <class KidsProc, class RawObs>
245void Beammap::pipeline(KidsProc &kidsproc, RawObs &rawobs) {
246 // only get kids params if not simulation
247 if (!telescope.sim_obs) {
248 // add kids models to apt
249 auto [kids_models, kids_model_header] = kidsproc.load_fit_report(rawobs);
250
251 Eigen::Index i = 0;
252 // loop through kids header
253 for (const auto &h: kids_model_header) {
254 std::string name = h;
255 if (name=="flag") {
256 name = "kids_flag";
257 }
258 calib.apt[name].resize(calib.n_dets);
259 Eigen::Index j = 0;
260 for (const auto &v: kids_models) {
261 calib.apt[name].segment(j,v.rows()) = v.col(i);
262 j = j + v.rows();
263 }
264
265 // search for key
266 bool found = false;
267 for (const auto &key: calib.apt_header_keys){
268 if (key==name) {
269 found = true;
270 }
271 }
272 // if not found, push back placeholder
273 if (!found) {
274 calib.apt_header_keys.push_back(name);
275 calib.apt_header_units[name] = "N/A";
276 }
277
278 // detector orientation
279 calib.apt_meta[name].push_back("units: N/A");
280 calib.apt_meta[name].push_back(name);
281 i++;
282 }
283 }
284
285 // run timestream pipeline
286 timestream_pipeline(kidsproc, rawobs);
287
288 // placeholder vectors of size nscans for grppi maps
289 scan_in_vec.resize(ptcs0.size());
290 std::iota(scan_in_vec.begin(), scan_in_vec.end(), 0);
291 scan_out_vec.resize(ptcs0.size());
292
293 // placeholder vectors of size ndet for grppi maps
294 det_in_vec.resize(n_maps);
295 std::iota(det_in_vec.begin(), det_in_vec.end(), 0);
296 det_out_vec.resize(n_maps);
297
298 // run iterative pipeline
300}
301
302template <class KidsProc, class RawObs>
303void Beammap::timestream_pipeline(KidsProc &kidsproc, RawObs &rawobs) {
304 using tuple_t = std::tuple<TCData<TCDataKind::RTC, Eigen::MatrixXd>,
305 std::vector<kids::KidsData<kids::KidsDataKind::RawTimeStream>>>;
306 // initialize number of completed scans
307 n_scans_done = 0;
308
309 // progress bar
310 tula::logging::progressbar pb(
311 [&](const auto &msg) { logger->info("{}", msg); }, 100, "RTC progress ");
312
313 // grppi generator function. gets time chunk data from files sequentially and passes them to grppi::farm
314 grppi::pipeline(tula::grppi_utils::dyn_ex(parallel_policy),
315 [&]() -> std::optional<tuple_t> {
316
317 // variable to hold current scan
318 static int scan = 0;
319 // loop through scans
320 while (scan < telescope.scan_indices.cols()) {
321 // update progress bar
322 pb.count(telescope.scan_indices.cols(), 1);
323
324 // create rtcdata
325 TCData<TCDataKind::RTC, Eigen::MatrixXd> rtcdata;
326 // get scan indices
327 rtcdata.scan_indices.data = telescope.scan_indices.col(scan);
328 // current scan
329 rtcdata.index.data = scan;
330
331 // vector to store kids data
332 std::vector<kids::KidsData<kids::KidsDataKind::RawTimeStream>> scan_rawobs;
333 // get kids data
334 if (!interp_over_gaps) {
335 scan_rawobs = kidsproc.load_rawobs(rawobs, scan, telescope.scan_indices, start_indices, end_indices);
336 }
337 else {
338 scan_rawobs = kidsproc.load_rawobs_gaps(rawobs, scan, telescope.scan_indices, start_indices,
339 t_common, nw_times, 1 / (2 * telescope.fsmp));
340 }
341 // current length of outer scans
342 Eigen::Index sl = rtcdata.scan_indices.data(3) - rtcdata.scan_indices.data(2) + 1;
343
344 // get raw tod from files
345 if (!interp_over_gaps) {
346 rtcdata.scans.data = kidsproc.populate_rtc(scan_rawobs, sl, calib.n_dets, tod_type);
347 }
348 else {
349 rtcdata.scans.data = kidsproc.populate_rtc_gaps(scan_rawobs, t_common, nw_times, masks, scan, 1 / (2 * telescope.fsmp),
351 }
352
353 // increment scan
354 scan++;
355 // return rtcdata, kidsproc, and raw data
356 return tuple_t(rtcdata,scan_rawobs);
357 }
358 // reset scan to zero for each obs
359 scan = 0;
360 return {};
361 },
362 // run the raw time chunk processing
363 run_timestream(kidsproc));
364}
365
366template <class KidsProc>
367auto Beammap::run_timestream(KidsProc &kidsproc) {
368 auto farm = grppi::farm(n_threads,[&](auto &input_tuple) -> TCData<TCDataKind::PTC,Eigen::MatrixXd> {
369 // RTCData input
370 auto& rtcdata = std::get<0>(input_tuple);
371 // start index input
372 auto& scan_rawobs = std::get<1>(input_tuple);
373
374 // allocate up bitwise timestream flags
375 rtcdata.flags2.data.setConstant(timestream::TimestreamFlags::Good);
376
377 // starting index for scan (outer scan)
378 Eigen::Index si = rtcdata.scan_indices.data(2);
379
380 // current length of outer scans
381 Eigen::Index sl = rtcdata.scan_indices.data(3) - rtcdata.scan_indices.data(2) + 1;
382
383 // copy scan's telescope vectors
384 for (const auto& x: telescope.tel_data) {
385 rtcdata.tel_data.data[x.first] = telescope.tel_data[x.first].segment(si,sl);
386 }
387
388 // copy pointing offsets
389 for (const auto& [axis,offset]: pointing_offsets_arcsec) {
390 rtcdata.pointing_offsets_arcsec.data[axis] = offset.segment(si,sl);
391 }
392
393 // get hwpr
395 rtcdata.hwpr_angle.data = calib.hwpr_angle.segment(si + hwpr_start_indices, sl);
396 }
397
398 // set up flags
399 rtcdata.flags.data.resize(rtcdata.scans.data.rows(), rtcdata.scans.data.cols());
400 rtcdata.flags.data.setConstant(0);
401
402 if (interp_over_gaps) {
403 int i = 0;
404 for (auto const& [key, val] : calib.nw_limits) {
405 auto& mask = masks[i];
406
407 Eigen::Index start = std::get<0>(calib.nw_limits[key]);
408 Eigen::Index end = std::get<1>(calib.nw_limits[key]) - 1;
409
410 for (int j = 0; j < rtcdata.flags.data.rows(); ++j) {
411 int start_index = j;
412 int size = 1;
414 start_index = std::max(j, static_cast<int>(j - rtcproc.filter.n_terms));
415 int end_index = std::min(j + rtcproc.filter.n_terms, rtcdata.flags.data.rows() - 1);
416 size = end_index - start_index + 1;
417 }
418 if (mask(j + si) == 0) {
419 rtcdata.flags.data.block(start_index, start, size, end - start + 1).setOnes();
420 }
421 }
422 logger->debug("{}/{} gaps flagged", rtcdata.flags.data.col(start).template cast<int>().sum(), rtcdata.flags.data.rows());
423 i++;
424 }
425 }
426
427 // create PTCData
429
430 logger->info("starting scan {}. {}/{} scans completed", rtcdata.index.data + 1, n_scans_done,
431 telescope.scan_indices.cols());
432
433 // run rtcproc
434 logger->info("raw time chunk processing for scan {}", rtcdata.index.data + 1);
435 auto map_indices = rtcproc.run(rtcdata, ptcdata, calib, telescope, omb.pixel_size_rad, map_grouping);
436
437 if (map_grouping!="detector") {
438 // remove flagged detectors
440 }
441
442 // remove outliers before cleaning
443 auto calib_scan = rtcproc.remove_bad_dets(ptcdata, calib, map_grouping);
444
445 // remove duplicate tones
446 if (!telescope.sim_obs) {
447 calib_scan = rtcproc.remove_nearby_tones(ptcdata, calib_scan, map_grouping);
448 }
449
450 // write rtc timestreams
451 if (run_tod_output && !tod_filename.empty()) {
452 if (tod_output_type == "rtc" || tod_output_type == "both") {
453 logger->info("writing raw time chunk");
455 ptcdata.pointing_offsets_arcsec.data, calib_scan);
456 }
457 }
458
459 // store indices for each ptcdata
460 ptcdata.map_indices.data = std::move(map_indices);
461
462 // move out ptcdata the PTCData vector at corresponding index
463 ptcs0.at(ptcdata.index.data) = std::move(ptcdata);
464 calib_scans0.at(ptcdata.index.data) = std::move(calib_scan);
465
466 // increment number of completed scans
467 n_scans_done++;
468 logger->info("done with scan {}. {}/{} scans completed", ptcdata.index.data + 1, n_scans_done, telescope.scan_indices.cols());
469
470 return ptcdata;
471 });
472
473 return farm;
474}
475
477 // run iterative stage
478 run_loop();
479
480 // write map summary
481 if (verbose_mode) {
483 }
484
485 // empty initial ptcdata vector to save memory
486 ptcs0.clear();
487
488 // set to input parallel policy
490
491 if (map_grouping=="detector") {
492 logger->info("calculating sensitivity");
493 // parallelize on detectors
494 grppi::map(tula::grppi_utils::dyn_ex(parallel_policy), det_in_vec, det_out_vec, [&](auto i) {
495 Eigen::MatrixXd det_sens, noise_flux;
496 // calc sensitivity within psd freq range
497 calc_sensitivity(ptcs, det_sens, noise_flux, telescope.d_fsmp, i, {sens_psd_limits_Hz(0), sens_psd_limits_Hz(1)});
498 // copy into apt table
499 calib.apt["sens"](i) = tula::alg::median(det_sens);
500
501 return 0;
502 });
503 }
504
505 // apt and sensitivity only relevant if beammapping
506 if (map_grouping=="detector") {
507 // rescale fit params from pixel to on-sky units
508 calib.apt["amp"] = params.col(0);
509 calib.apt["x_t"] = RAD_TO_ASEC*omb.pixel_size_rad*(params.col(1).array() - (omb.n_cols)/2);
510 calib.apt["y_t"] = RAD_TO_ASEC*omb.pixel_size_rad*(params.col(2).array() - (omb.n_rows)/2);
513 calib.apt["angle"] = params.col(5);
514
515 // rescale fit errors from pixel to on-sky units
516 calib.apt["amp_err"] = perrors.col(0);
517 calib.apt["x_t_err"] = RAD_TO_ASEC*omb.pixel_size_rad*(perrors.col(1));
518 calib.apt["y_t_err"] = RAD_TO_ASEC*omb.pixel_size_rad*(perrors.col(2));
519 calib.apt["a_fwhm_err"] = RAD_TO_ASEC*STD_TO_FWHM*omb.pixel_size_rad*(perrors.col(3));
520 calib.apt["b_fwhm_err"] = RAD_TO_ASEC*STD_TO_FWHM*omb.pixel_size_rad*(perrors.col(4));
521 calib.apt["angle_err"] = perrors.col(5);
522
523 // add convergence iteration to apt table
524 calib.apt["converge_iter"] = converge_iter.cast<double> ();
525
526 // flag detectors in apt based on config limits
528
529 // subtract reference detector position and derotate
530 process_apt();
531
532 // add final apt table to timestream files
533 if (run_tod_output && !tod_filename.empty()) {
534 // vectors to hold tangent plane pointing for all ptcs (n_chunks x [n_pts x n_dets])
535 std::vector<Eigen::MatrixXd> lat, lon;
536
537 // recalculate tangent plane pointing for tod output
538 for (Eigen::Index i=0; i<ptcs.size(); ++i) {
539 // tangent plane pointing for each detector
540 Eigen::MatrixXd ptc_lat(ptcs[i].scans.data.rows(), ptcs[i].scans.data.cols());
541 Eigen::MatrixXd ptc_lon(ptcs[i].scans.data.rows(), ptcs[i].scans.data.cols());
542 // loop through detectors
543 grppi::map(tula::grppi_utils::dyn_ex(parallel_policy), det_in_vec, det_out_vec, [&](auto j) {
544 // det indices
545 auto det_index = j;
546 double az_off = calib.apt["x_t"](det_index);
547 double el_off = calib.apt["y_t"](det_index);
548
549 // get tangent pointing
550 auto [det_lat, det_lon] = engine_utils::calc_det_pointing(ptcs[i].tel_data.data, az_off,
551 el_off, telescope.pixel_axes,
552 ptcs[i].pointing_offsets_arcsec.data,
554 ptc_lat.col(j) = std::move(det_lat);
555 ptc_lon.col(j) = std::move(det_lon);
556
557 return 0;
558 });
559 lat.push_back(std::move(ptc_lat));
560 lon.push_back(std::move(ptc_lon));
561 }
562
563 // set parallel policy to sequential for tod output
564 parallel_policy = "seq";
565
566 logger->info("adding final apt and detector pointing to tod files");
567 // loop through tod files
568 for (const auto & [key, val]: tod_filename) {
569 netCDF::NcFile fo(val, netCDF::NcFile::write);
570 // overwrite apt table
571 for (auto const& x: calib.apt) {
572 if (x.first!="flag2") {
573 // start index for apt table
574 std::vector<std::size_t> start_index_apt = {0};
575 // size for apt
576 std::vector<std::size_t> size_apt = {1};
577 netCDF::NcVar apt_v = fo.getVar("apt_" + x.first);
578 if (!apt_v.isNull()) {
579 for (std::size_t i=0; i< TULA_SIZET(calib.n_dets); ++i) {
580 start_index_apt[0] = i;
581 apt_v.putVar(start_index_apt, size_apt, &calib.apt[x.first](i));
582 }
583 }
584 }
585 }
586
587 // detector tangent plane pointing
588 netCDF::NcVar det_lat_v = fo.getVar("det_lat");
589 netCDF::NcVar det_lon_v = fo.getVar("det_lon");
590
591 // detector absolute pointing
592 netCDF::NcVar det_ra_v = fo.getVar("det_ra");
593 netCDF::NcVar det_dec_v = fo.getVar("det_dec");
594
595 // start indices for data
596 std::vector<std::size_t> start_index = {0, 0};
597 // size for data
598 std::vector<std::size_t> size = {1, TULA_SIZET(calib.n_dets)};
599 std::size_t k = 0;
600 // loop through ptcs
601 for (Eigen::Index i=0; i<lat.size(); ++i) {
602 // loop through n_pts
603 for (std::size_t j=0; j < TULA_SIZET(lat[i].rows()); ++j) {
604 start_index[0] = k;
605 k++;
606 // append detector latitudes
607 Eigen::VectorXd lat_row = lat[i].row(j);
608 det_lat_v.putVar(start_index, size, lat_row.data());
609
610 // append detector longitudes
611 Eigen::VectorXd lon_row = lon[i].row(j);
612 det_lon_v.putVar(start_index, size, lon_row.data());
613
614 if (telescope.pixel_axes == "radec") {
615 // get absolute pointing
616 auto [dec, ra] = engine_utils::tangent_to_abs(lat_row, lon_row, telescope.tel_header["Header.Source.Ra"](0),
617 telescope.tel_header["Header.Source.Dec"](0));
618 // append detector ra
619 det_ra_v.putVar(start_index, size, ra.data());
620
621 // append detector dec
622 det_dec_v.putVar(start_index, size, dec.data());
623 }
624 }
625 }
626 }
627
628 // empty ptcdata vector to save memory
629 ptcs.clear();
630 }
631 }
632
633 else {
634 // calculate map psds
635 logger->info("calculating map psd");
637 // calculate map histograms
638 logger->info("calculating map histogram");
640 }
641}
642
643
645 // variable to control iteration
646 bool keep_going = true;
647
648 // declare random number generator
649 boost::random::mt19937 eng;
650
651 // boost random number generator (0,1)
652 boost::random::uniform_int_distribution<> rands{0,1};
653
654 // iterative loop
655 while (keep_going) {
656 logger->info("starting iter {}", current_iter);
657
658 // copy ptcs
659 ptcs = ptcs0;
660 // copy calibs
662
663 // copy signal for convergence test
666 // calc mean rms
667 if (current_iter == 1) {
668 // use obs map buffer
669 if (!omb.noise.empty()) {
671 }
672 }
673 }
674
675 // progress bar
676 tula::logging::progressbar pb(
677 [&](const auto &msg) { logger->info("{}", msg); }, 100, "PTC progress ");
678
679
680 // cleaning (separate from mapmaking loop due to jinc mapmaking parallelization)
681 grppi::map(tula::grppi_utils::dyn_ex(omb.parallel_policy), scan_in_vec, scan_out_vec, [&](auto i) {
682 if (run_mapmaking) {
683 if (current_iter > 0) {
684 if (!ptcproc.run_fruit_loops) {
685 // if not running fruit loops use source fit
686 logger->info("subtracting gaussian from tod");
687 // subtract gaussian
688 ptcproc.add_gaussian<timestream::TCProc::SourceType::NegativeGaussian>(ptcs[i], params, telescope.pixel_axes, map_grouping,
689 calib.apt,omb.pixel_size_rad, omb.n_rows, omb.n_cols);
690 }
691 else {
692 logger->info("subtracting map from tod");
693 // subtract map
694 ptcproc.map_to_tod<timestream::TCProc::SourceType::NegativeMap>(omb, ptcs[i], calib, ptcs[i].map_indices.data, telescope.pixel_axes,
695 map_grouping);
696 }
697 }
698 }
699
700 // clean the maps
701 logger->info("processed time chunk processing for scan {}", i + 1);
703
704 if (run_mapmaking) {
705 if (current_iter > 0) {
706 // if not running fruit loops use source fit
708 logger->info("adding gaussian to tod");
709 // add gaussian back
712 }
713 else {
714 logger->info("adding map to tod");
715 // add map back
718 }
719 }
720 }
721
722 // remove outliers after clean
724
725 if (map_grouping == "detector") {
726 // set weights to a constant value
727 ptcs[i].weights.data.resize(ptcs[i].scans.data.cols());
728 ptcs[i].weights.data.setOnes();
729 }
730 else {
731 // calculate weights
732 logger->info("calculating weights");
734
735 // reset weights to median
737 }
738
739 // write out chunk summary
741 logger->debug("writing chunk summary");
743 }
744
745 // calc stats
746 logger->debug("calculating stats");
748
749 return 0;
750 });
751
752 // write ptc timestreams
753 if (run_tod_output && !tod_filename.empty()) {
754 if (tod_output_type == "ptc" || tod_output_type == "both") {
755 logger->info("writing processed time chunk");
757 for (Eigen::Index i=0; i<telescope.scan_indices.cols(); ++i) {
759 ptcs[i].pointing_offsets_arcsec.data, calib_scans[i]);
760 }
761 }
762 }
763 }
764
765 logger->info("starting mapmaking");
766
767 if (run_mapmaking) {
768 // set maps to zero for each iteration
769 for (Eigen::Index i=0; i<n_maps; ++i) {
770 omb.signal[i].setZero();
771 omb.weight[i].setZero();
772
773 // clear coverage
774 if (!omb.coverage.empty()) {
775 omb.coverage[i].setZero();
776 }
777 // clear kernel
778 if (rtcproc.run_kernel) {
779 omb.kernel[i].setZero();
780 }
781 // clear noise
782 if (!omb.noise.empty()) {
783 omb.noise[i].setZero();
784 }
785
786 if (run_noise) {
787 for (auto& ptcdata: ptcs) {
788 if (omb.randomize_dets) {
789 ptcdata.noise.data = Eigen::Matrix<int, Eigen::Dynamic, Eigen::Dynamic>::Zero(omb.n_noise, calib.n_dets)
790 .unaryExpr([&](int dummy){ return 2 * rands(eng) - 1; });
791 } else {
792 ptcdata.noise.data = Eigen::Matrix<int, Eigen::Dynamic, 1>::Zero(omb.n_noise)
793 .unaryExpr([&](int dummy){ return 2 * rands(eng) - 1; });
794 }
795 }
796 }
797 }
798
799 logger->info("running mapmaking");
800
801 if (map_grouping == "detector") {
802 bool run_omb = true;
803 for (auto& ptc : ptcs) {
804 if (map_method == "naive") {
805 // naive mapmaker
807 calib.apt, telescope.d_fsmp, run_omb, run_noise);
808 }
809 else if (map_method == "jinc") {
810 // jinc mapmaker
811 jinc_mm.populate_maps_jinc_parallel(ptc, omb, cmb, ptc.map_indices.data, telescope.pixel_axes,
812 calib.apt, telescope.d_fsmp, run_omb, run_noise);
813 }
814 // update progress bar
815 pb.count(telescope.scan_indices.cols(), 1);
816 }
817 }
818
819 else {
820 // mapmaking
821 grppi::map(tula::grppi_utils::dyn_ex(map_parallel_policy), scan_in_vec, scan_out_vec, [&](auto i) {
822 bool run_omb = true;
823 // naive mapmaker
824 if (map_method == "naive") {
825 naive_mm.populate_maps_naive(ptcs[i], omb, cmb, ptcs[i].map_indices.data, telescope.pixel_axes,
826 calib.apt, telescope.d_fsmp, run_omb, run_noise);
827 }
828 else if (map_method == "jinc") {
829 jinc_mm.populate_maps_jinc(ptcs[i], omb, cmb, ptcs[i].map_indices.data, telescope.pixel_axes,
830 calib.apt, telescope.d_fsmp, run_omb, run_noise);
831 }
832
833 // update progress bar
834 pb.count(telescope.scan_indices.cols(), 1);
835
836 return 0;
837 });
838 }
839
840 // normalize maps
841 logger->info("normalizing maps");
843
844 // initial position for fitting
845 double init_row = -99;
846 double init_col = -99;
847
848 logger->info("fitting maps");
849 grppi::map(tula::grppi_utils::dyn_ex(omb.parallel_policy), det_in_vec, det_out_vec, [&](auto i) {
850 // only fit if not converged
851 if (!converged(i)) {
852 // get array number
853 auto array = maps_to_arrays(i);
854 // get initial guess fwhm from theoretical fwhms for the arrays
855 double init_fwhm = toltec_io.array_fwhm_arcsec[array]*ASEC_TO_RAD/omb.pixel_size_rad;
856 // fit the maps
857 auto [det_params, det_perror, good_fit] =
858 map_fitter.fit_to_gaussian<engine_utils::mapFitter::beammap>(omb.signal[i], omb.weight[i],
859 init_fwhm, init_row, init_col);
860
861 params.row(i) = det_params;
862 perrors.row(i) = det_perror;
863 good_fits(i) = good_fit;
864 }
865 // otherwise keep value from previous iteration
866 else {
867 params.row(i) = p0.row(i);
868 perrors.row(i) = perror0.row(i);
869 }
870
871 return 0;
872 });
873
874 logger->info("number of good fits {}/{}", good_fits.cast<double>().sum(), n_maps);
875 }
876
877 // increment loop iteration
878 current_iter++;
879
880 if (current_iter < beammap_iter_max) {
881 // check if all detectors are converged
882 if ((converged.array() == true).all()) {
883 logger->info("all maps converged");
884 keep_going = false;
885 }
886 else if (current_iter > 1) {
887 // only do convergence test if tolerance is above zero, otherwise run all iterations
888 if (beammap_iter_tolerance > 0) {
889 // loop through maps and check if it is converged
890 logger->info("checking convergence");
891 grppi::map(tula::grppi_utils::dyn_ex(omb.parallel_policy), det_in_vec, det_out_vec, [&](auto i) {
892 if (!converged(i)) {
893 // get relative change from last iteration
894 Eigen::ArrayXd diff;
895 if (!ptcproc.run_fruit_loops) {
896 diff = abs((params.row(i).array() - p0.row(i).array())/p0.row(i).array());
897 }
898 else {
899 diff = abs(omb_copy.signal[i].array() - omb.signal[i].array()/omb_copy.signal[i].array());
900 }
901 // if a variable is constant, make sure no nans are present
902 auto d = (diff.array()).isNaN().select(0,diff);
903 if ((d.array() <= beammap_iter_tolerance).all()) {
904 // set as converged
905 converged(i) = true;
906 // set convergence iteration
907 converge_iter(i) = current_iter;
908 }
909 }
910 return 0;
911 });
912
913 logger->info("{} maps converged on iter {}", (converged.array() == true).count(), current_iter);
914
915 // stop if all maps converged
916 if ((converged.array() == true).all()) {
917 logger->info("all maps converged");
918 keep_going = false;
919 }
920 }
921 else {
922 logger->info("bypassing convergence check");
923 }
924 }
925
926 // set previous iteration fits to current iteration fits
927 p0 = params;
928 perror0 = perrors;
929 }
930 else {
931 logger->info("max iteration reached");
932 keep_going = false;
933 }
934 }
935}
936
938 // setup bitwise flags
939 flag2.resize(calib.n_dets);
940 flag2.setConstant(AptFlags::Good);
941
942 // track number of flagged detectors
943 int n_flagged_dets = 0;
944
945 logger->info("flagging detectors");
946 // first flag based on fit values and signal-to-noise
947 grppi::map(tula::grppi_utils::dyn_ex(parallel_policy), det_in_vec, det_out_vec, [&](auto i) {
948 // get array of current detector
949 auto array_index = calib.apt["array"](i);
950 std::string array_name = toltec_io.array_name_map[array_index];
951
952 // calculate map standard deviation
953 double map_std_dev = engine_utils::calc_std_dev(omb.signal[i]);
954
955 // set apt signal to noise
956 calib.apt["sig2noise"](i) = params(i,0)/perrors(i,0);
957
958 // flag bad fits
959 if (!good_fits(i)) {
960 if (calib.apt["flag"](i)==0) {
961 n_flagged_dets++;
962 calib.apt["flag"](i) = 1;
963 }
965 }
966 // flag detectors with outlier a_fwhm values
967 if (calib.apt["a_fwhm"](i) < lower_fwhm_arcsec[array_name] ||
968 ((calib.apt["a_fwhm"](i) > upper_fwhm_arcsec[array_name]) && upper_fwhm_arcsec[array_name] > 0)) {
969 if (calib.apt["flag"](i)==0) {
970 n_flagged_dets++;
971 calib.apt["flag"](i) = 1;
972 }
974 }
975 // flag detectors with outlier b_fwhm values
976 if (calib.apt["b_fwhm"](i) < lower_fwhm_arcsec[array_name] ||
977 ((calib.apt["b_fwhm"](i) > upper_fwhm_arcsec[array_name] && upper_fwhm_arcsec[array_name] > 0))) {
978 if (calib.apt["flag"](i)==0) {
979 n_flagged_dets++;
980 calib.apt["flag"](i) = 1;
981 }
983 }
984 // flag detectors with outlier S/N values
985 if ((params(i,0)/map_std_dev < lower_sig2noise[array_name]) ||
986 ((params(i,0)/map_std_dev > upper_sig2noise[array_name]) && (upper_sig2noise[array_name] > 0))) {
987 if (calib.apt["flag"](i)==0) {
988 n_flagged_dets++;
989 calib.apt["flag"](i) = 1;
990 }
992 }
993 return 0;
994 });
995
996
997 // median network sensitivity for flagging
998 std::map<Eigen::Index, double> nw_median_sens;
999
1000 // calc median sens from unflagged detectors for each nw
1001 logger->debug("calculating mean sensitivities");
1002 for (Eigen::Index i=0; i<calib.n_nws; ++i) {
1003 Eigen::Index nw = calib.nws(i);
1004
1005 // nw sensitivity
1006 auto nw_sens = calib.apt["sens"](Eigen::seq(std::get<0>(calib.nw_limits[nw]),
1007 std::get<1>(calib.nw_limits[nw])-1));
1008
1009 // number of good detectors
1010 Eigen::Index n_good_det = (calib.apt["flag"](Eigen::seq(std::get<0>(calib.nw_limits[nw]),
1011 std::get<1>(calib.nw_limits[nw])-1)).array()==0).count();
1012
1013 if (n_good_det>0) {
1014 // to hold good detectors
1015 Eigen::VectorXd sens(n_good_det);
1016
1017 // remove flagged dets
1018 Eigen::Index j = std::get<0>(calib.nw_limits[nw]);
1019 Eigen::Index k = 0;
1020 for (Eigen::Index m=0; m<sens.size(); m++) {
1021 if (calib.apt["flag"](j)==0) {
1022 sens(k) = nw_sens(m);
1023 k++;
1024 }
1025 j++;
1026 }
1027 // calculate median sens
1028 nw_median_sens[nw] = tula::alg::median(sens);
1029 }
1030 else {
1031 nw_median_sens[nw] = tula::alg::median(nw_sens);
1032 }
1033 }
1034
1035
1036 // flag too low/high sensitivies based on the median unflagged sensitivity of each nw
1037 logger->debug("flagging sensitivities");
1038 grppi::map(tula::grppi_utils::dyn_ex(parallel_policy), det_in_vec, det_out_vec, [&](auto i) {
1039 // get nw of current detector
1040 auto nw_index = calib.apt["nw"](i);
1041
1042 // flag outlier sensitivities
1043 if (calib.apt["sens"](i) < lower_sens_factor*nw_median_sens[nw_index] ||
1044 (calib.apt["sens"](i) > upper_sens_factor*nw_median_sens[nw_index] && upper_sens_factor > 0)) {
1045 if (calib.apt["flag"](i)==0) {
1046 calib.apt["flag"](i) = 1;
1047 n_flagged_dets++;
1048 }
1049 flag2(i) |= AptFlags::Sens;
1050 }
1051
1052 return 0;
1053 });
1054
1055 // std maps to hold median unflagged x and y positions
1056 std::map<std::string, double> array_median_x_t, array_median_y_t;
1057
1058 // calc median x_t and y_t values from unflagged detectors for each arrays
1059 logger->debug("calculating array median positions");
1060 for (Eigen::Index i=0; i<calib.n_arrays; ++i) {
1061 Eigen::Index array = calib.arrays(i);
1062 std::string array_name = toltec_io.array_name_map[array];
1063
1064 // x_t
1065 auto array_x_t = calib.apt["x_t"](Eigen::seq(std::get<0>(calib.array_limits[array]),
1066 std::get<1>(calib.array_limits[array])-1));
1067 // y_t
1068 auto array_y_t = calib.apt["y_t"](Eigen::seq(std::get<0>(calib.array_limits[array]),
1069 std::get<1>(calib.array_limits[array])-1));
1070 // number of good detectors
1071 Eigen::Index n_good_det = (calib.apt["flag"](Eigen::seq(std::get<0>(calib.array_limits[array]),
1072 std::get<1>(calib.array_limits[array])-1)).array()==0).count();
1073
1074 // to hold good detectors
1075 Eigen::VectorXd x_t, y_t;
1076
1077 if (n_good_det>0) {
1078 x_t.resize(n_good_det);
1079 y_t.resize(n_good_det);
1080
1081 // remove flagged dets
1082 Eigen::Index j = std::get<0>(calib.array_limits[array]);
1083 Eigen::Index k = 0;
1084 for (Eigen::Index m=0; m<array_x_t.size(); m++) {
1085 if (calib.apt["flag"](j)==0) {
1086 x_t(k) = array_x_t(m);
1087 y_t(k) = array_y_t(m);
1088 k++;
1089 }
1090 j++;
1091 }
1092 // calculate medians
1093 array_median_x_t[array_name] = tula::alg::median(x_t);
1094 array_median_y_t[array_name] = tula::alg::median(y_t);
1095 }
1096 else {
1097 // if no good dets, use all dets to calculate median
1098 array_median_x_t[array_name] = tula::alg::median(array_x_t);
1099 array_median_y_t[array_name] = tula::alg::median(array_y_t);
1100 }
1101 }
1102
1103 // remove detectors above distance limits
1104 logger->debug("flagging detector positions");
1105 grppi::map(tula::grppi_utils::dyn_ex(parallel_policy), det_in_vec, det_out_vec, [&](auto i) {
1106 // get array of current detector
1107 auto array_index = calib.apt["array"](i);
1108 std::string array_name = toltec_io.array_name_map[array_index];
1109
1110 // calculate distance of detector from mean position of all detectors
1111 double dist = sqrt(pow(calib.apt["x_t"](i) - array_median_x_t[array_name],2) +
1112 pow(calib.apt["y_t"](i) - array_median_y_t[array_name],2));
1113
1114 // flag detectors that are further than the mean value than the distance limit
1115 if (dist > max_dist_arcsec[array_name] && max_dist_arcsec[array_name] > 0) {
1116 if (calib.apt["flag"](i)==0) {
1117 n_flagged_dets++;
1118 calib.apt["flag"](i) = 1;
1119 }
1121 }
1122
1123 return 0;
1124 });
1125
1126 // print number of flagged detectors
1127 logger->info("{} detectors were flagged", n_flagged_dets);
1128
1129 // calculate fcf
1130 logger->debug("calculating flux conversion factors");
1131 grppi::map(tula::grppi_utils::dyn_ex(parallel_policy), det_in_vec, det_out_vec, [&](auto i) {
1132 // get array of current detector
1133 auto array_index = calib.apt["array"](i);
1134 std::string array_name = toltec_io.array_name_map[array_index];
1135
1136 // calc flux scale (always in mJy/beam)
1137 if (params(i,0)!=0) {
1138 calib.apt["flxscale"](i) = beammap_fluxes_mJy_beam[array_name]/params(i,0);
1139 calib.apt["sens"](i) = calib.apt["sens"](i)*calib.apt["flxscale"](i);
1140 }
1141 // set fluxscale (fcf) to zero if flagged
1142 else {
1143 calib.apt["flxscale"](i) = 0;
1144 calib.apt["sens"](i) = 0;
1145 }
1146 return 0;
1147 });
1148
1149 // re-run calib setup to get average fwhms and beam areas
1150 calib.setup();
1151
1152 // calculate source flux in MJy/sr from average beamsizes
1153 for (Eigen::Index i=0; i<calib.n_arrays; ++i) {
1154 Eigen::Index array = calib.arrays(i);
1155 std::string array_name = toltec_io.array_name_map[array];
1156
1157 // get source flux in MJy/Sr
1159 }
1160}
1161
1163 // reference detector x and y
1164 double ref_det_x_t = 0;
1165 double ref_det_y_t = 0;
1166
1167 // initial reference det
1169
1170 // if particular reference detector is requested
1172 if (beammap_reference_det >= 0) {
1174 // set reference x_t and y_t
1175 ref_det_x_t = calib.apt["x_t"](beammap_reference_det_found);
1176 ref_det_y_t = calib.apt["y_t"](beammap_reference_det_found);
1177 }
1178 // else use closest to (0,0) in map 0 apt
1179 else {
1180 logger->info("finding a reference detector");
1181 // calc x_t and y_t values from unflagged detectors for each arrays
1182 Eigen::Index array = calib.arrays(0);
1183
1184 // x_t
1185 auto array_x_t = calib.apt["x_t"](Eigen::seq(std::get<0>(calib.array_limits[array]),
1186 std::get<1>(calib.array_limits[array])-1));
1187 // y_t
1188 auto array_y_t = calib.apt["y_t"](Eigen::seq(std::get<0>(calib.array_limits[array]),
1189 std::get<1>(calib.array_limits[array])-1));
1190 // number of good detectors
1191 Eigen::Index n_good_det = (calib.apt["flag"](Eigen::seq(std::get<0>(calib.array_limits[array]),
1192 std::get<1>(calib.array_limits[array])-1)).array()==0).count();
1193
1194 Eigen::VectorXd x_t, y_t, det_indices, dist;
1195
1196 if (n_good_det>0) {
1197 x_t.resize(n_good_det);
1198 y_t.resize(n_good_det);
1199 det_indices.resize(n_good_det);
1200
1201 // get good detector positions
1202 Eigen::Index j = std::get<0>(calib.array_limits[array]);
1203 Eigen::Index k = 0;
1204 for (Eigen::Index i=0; i<array_x_t.size(); ++i) {
1205 if (calib.apt["flag"](j)==0) {
1206 x_t(k) = array_x_t(i);
1207 y_t(k) = array_y_t(i);
1208 det_indices(k) = j;
1209 k++;
1210 }
1211 j++;
1212 }
1213
1214 // distance from (0,0)
1215 dist = pow(x_t.array(),2) + pow(y_t.array(),2);
1216 }
1217 else {
1218 dist = pow(array_x_t.array(),2) + pow(array_y_t.array(),2);
1219 det_indices = Eigen::VectorXd::LinSpaced(array_x_t.size(), std::get<0>(calib.array_limits[array]),
1220 std::get<1>(calib.array_limits[array]));
1221 }
1222
1223 // index of detector closest to zero
1224 auto min_dist = dist.minCoeff(&beammap_reference_det_found);
1225
1226 // get row in apt table
1228
1229 // set reference x_t and y_t
1230 ref_det_x_t = calib.apt["x_t"](beammap_reference_det_found);
1231 ref_det_y_t = calib.apt["y_t"](beammap_reference_det_found);
1232 }
1233 logger->info("using detector {} at ({},{}) arcsec",beammap_reference_det_found,
1234 static_cast<float>(ref_det_x_t),static_cast<float>(ref_det_y_t));
1235 }
1236 else {
1237 logger->info("no reference detector selected");
1238 }
1239
1240 // add reference detector to APT meta data
1241 calib.apt_meta["reference_x_t"] = ref_det_x_t;
1242 calib.apt_meta["reference_y_t"] = ref_det_y_t;
1243
1244 // raw (not derotated or reference detector subtracted) detector x and y values
1245 calib.apt["x_t_raw"] = calib.apt["x_t"];
1246 calib.apt["y_t_raw"] = calib.apt["y_t"];
1247
1248 // align to reference detector if specified and subtract its position from x and y
1249 calib.apt["x_t"] = calib.apt["x_t"].array() - ref_det_x_t;
1250 calib.apt["y_t"] = calib.apt["y_t"].array() - ref_det_y_t;
1251
1252 // derotated detector x and y values
1253 calib.apt["x_t_derot"] = calib.apt["x_t"];
1254 calib.apt["y_t_derot"] = calib.apt["y_t"];
1255
1256 // derotation elevation
1257 calib.apt["derot_elev"].setConstant(telescope.tel_data["TelElAct"].mean());
1258
1259 // calculate derotated positions
1260 Eigen::VectorXd rot_az_off = cos(-calib.apt["derot_elev"].array())*calib.apt["x_t_derot"].array() -
1261 sin(-calib.apt["derot_elev"].array())*calib.apt["y_t_derot"].array();
1262 Eigen::VectorXd rot_alt_off = sin(-calib.apt["derot_elev"].array())*calib.apt["x_t_derot"].array() +
1263 cos(-calib.apt["derot_elev"].array())*calib.apt["y_t_derot"].array();
1264
1265 // overwrite x_t and y_t
1266 calib.apt["x_t_derot"] = -rot_az_off;
1267 calib.apt["y_t_derot"] = -rot_alt_off;
1268
1269 if (beammap_derotate) {
1270 logger->info("derotating apt");
1271 // if derotation requested set default positions to derotated positions
1272 calib.apt["x_t"] = calib.apt["x_t_derot"];
1273 calib.apt["y_t"] = calib.apt["y_t_derot"];
1274 }
1275}
1276
1277template <mapmaking::MapType map_type>
1279 // pointer to map buffer
1280 mapmaking::MapBuffer* mb = nullptr;
1281 // pointer to data file fits vector
1282 std::vector<fitsIO<file_type_enum::write_fits, CCfits::ExtHDU*>>* f_io = nullptr;
1283 // pointer to noise file fits vector
1284 std::vector<fitsIO<file_type_enum::write_fits, CCfits::ExtHDU*>>* n_io = nullptr;
1285
1286 // directory name
1287 std::string dir_name;
1288
1289 // raw obs maps
1290 if constexpr (map_type == mapmaking::RawObs) {
1291 mb = &omb;
1292 f_io = &fits_io_vec;
1293 n_io = &noise_fits_io_vec;
1294 dir_name = obsnum_dir_name + "raw/";
1295
1296 // write stats file
1297 write_stats();
1298
1299 // add header informqtion to tod
1300 if (run_tod_output && !tod_filename.empty()) {
1302 }
1303
1304 // only write apt table if beammapping
1305 if (map_grouping=="detector") {
1306 logger->info("writing apt table");
1310
1311 Eigen::MatrixXd apt_table(calib.n_dets, calib.apt_header_keys.size());
1312
1313 // copy to table
1314 Eigen::Index i = 0;
1315 for (auto const& x: calib.apt_header_keys) {
1316 if (x != "flag2") {
1317 apt_table.col(i) = calib.apt[x];
1318 }
1319 else {
1320 apt_table.col(i) = flag2.cast<double> ();
1321 }
1322 i++;
1323 }
1324
1325 // write to ecsv
1326 to_ecsv_from_matrix(apt_filename, apt_table, calib.apt_header_keys, calib.apt_meta);
1327
1328 logger->info("done writing apt table {}.ecsv",apt_filename);
1329 }
1330 }
1331
1332 // filtered obs maps
1333 else if constexpr (map_type == mapmaking::FilteredObs) {
1334 mb = &omb;
1335 f_io = &filtered_fits_io_vec;
1337 dir_name = obsnum_dir_name + "filtered/";
1338 }
1339
1340 // raw coadded maps
1341 else if constexpr (map_type == mapmaking::RawCoadd) {
1342 mb = &cmb;
1343 f_io = &coadd_fits_io_vec;
1345 dir_name = coadd_dir_name + "raw/";
1346 }
1347
1348 // filtered coadded maps
1349 else if constexpr (map_type == mapmaking::FilteredCoadd) {
1350 mb = &cmb;
1353 dir_name = coadd_dir_name + "filtered/";
1354 }
1355
1356 if (run_mapmaking) {
1357 // wiener filtered maps write before this and are deleted from the vector.
1358 if (!f_io->empty()) {
1359 {
1360 // progress bar
1361 tula::logging::progressbar pb(
1362 [&](const auto &msg) { logger->info("{}", msg); }, 100, "output progress ");
1363
1364 for (Eigen::Index i=0; i<f_io->size(); ++i) {
1365 // get the array for the given map
1366 // add primary hdu
1367 logger->debug("adding primary header to file {}",i);
1368 add_phdu(f_io, mb, i);
1369
1370 if (!mb->noise.empty()) {
1371 logger->debug("adding primary header to noise file {}",i);
1372 add_phdu(n_io, mb, i);
1373 }
1374 }
1375
1376 logger->debug("done adding primary headers");
1377
1378 // write the maps
1379 Eigen::Index k = 0;
1380 Eigen::Index step = 2;
1381
1382 if (!mb->kernel.empty()) {
1383 step++;
1384 }
1385 if (!mb->coverage.empty()) {
1386 step++;
1387 }
1388
1389 // write the maps
1390 for (Eigen::Index i=0; i<n_maps; ++i) {
1391 // update progress bar
1392 pb.count(n_maps, 1);
1393 logger->debug("adding map");
1394 write_maps(f_io,n_io,mb,i);
1395
1396 if (map_grouping=="detector") {
1397 if constexpr (map_type == mapmaking::RawObs) {
1398 // get the array for the given map
1399 Eigen::Index map_index = arrays_to_maps(i);
1400
1401 // check if we move from one file to the next
1402 // if so go back to first hdu layer
1403 if (i>0) {
1404 if (map_index > arrays_to_maps(i-1)) {
1405 k = 0;
1406 }
1407 }
1408
1409 // add apt table
1410 logger->debug("adding beammap header keys");
1411 for (auto const& key: calib.apt_header_keys) {
1412 if (key!="flag2") {
1413 try {
1414 f_io->at(map_index).hdus.at(k)->addKey("BEAMMAP." + key, calib.apt[key](i), key
1415 + " (" + calib.apt_header_units[key] + ")");
1416 } catch(...) {
1417 f_io->at(map_index).hdus.at(k)->addKey("BEAMMAP." + key, 0.0, key
1418 + " (" + calib.apt_header_units[key] + ")");
1419 }
1420 }
1421 else {
1422 f_io->at(map_index).hdus.at(k)->addKey("BEAMMAP." + key, flag2(i), key
1423 + " (" + calib.apt_header_units[key] + ")");
1424 }
1425 }
1426 // increment hdu layer
1427 k = k + step;
1428 }
1429 }
1430 }
1431 }
1432
1433 logger->info("maps have been written to:");
1434 for (Eigen::Index i=0; i<f_io->size(); ++i) {
1435 logger->info("{}.fits",f_io->at(i).filepath);
1436 }
1437 }
1438
1439 // clear fits file vectors to ensure its closed.
1440 f_io->clear();
1441 n_io->clear();
1442
1443 if (map_grouping!="detector") {
1444 // write psd and histogram files
1445 logger->debug("writing psds");
1446 write_psd<map_type>(mb, dir_name);
1447 logger->debug("writing histograms");
1448 write_hist<map_type>(mb, dir_name);
1449 }
1450 }
1451}
1452
Definition beammap.h:16
Eigen::Vector< int, Eigen::Dynamic > converge_iter
Definition beammap.h:37
void run_loop()
Definition beammap.h:644
Eigen::Index beammap_reference_det_found
Definition beammap.h:46
Eigen::Matrix< uint16_t, Eigen::Dynamic, 1 > flag2
Definition beammap.h:60
void setup()
Definition beammap.h:101
Eigen::MatrixXd perror0
Definition beammap.h:40
AptFlags
Definition beammap.h:49
@ Position
Definition beammap.h:56
@ ElFWHM
Definition beammap.h:53
@ Good
Definition beammap.h:50
@ Sens
Definition beammap.h:55
@ AzFWHM
Definition beammap.h:52
@ Sig2Noise
Definition beammap.h:54
@ BadFit
Definition beammap.h:51
void timestream_pipeline(KidsProc &, RawObs &)
Definition beammap.h:303
void pipeline(KidsProc &, RawObs &)
Definition beammap.h:245
std::vector< int > det_in_vec
Definition beammap.h:67
std::vector< TCData< TCDataKind::PTC, Eigen::MatrixXd > > ptcs
Definition beammap.h:22
Eigen::Index current_iter
Definition beammap.h:31
void output()
Definition beammap.h:1278
std::vector< int > det_out_vec
Definition beammap.h:67
void loop_pipeline()
Definition beammap.h:476
Eigen::Matrix< bool, Eigen::Dynamic, 1 > good_fits
Definition beammap.h:63
std::vector< engine::Calib > calib_scans0
Definition beammap.h:28
std::vector< TCData< TCDataKind::PTC, Eigen::MatrixXd > > ptcs0
Definition beammap.h:22
Eigen::MatrixXd p0
Definition beammap.h:40
std::vector< int > scan_out_vec
Definition beammap.h:66
std::vector< engine::Calib > calib_scans
Definition beammap.h:28
Eigen::Matrix< bool, Eigen::Dynamic, 1 > converged
Definition beammap.h:34
Eigen::MatrixXd perrors
Definition beammap.h:43
std::string map_parallel_policy
Definition beammap.h:19
std::vector< int > scan_in_vec
Definition beammap.h:66
void set_apt_flags()
Definition beammap.h:937
mapmaking::MapBuffer omb_copy
Definition beammap.h:25
void process_apt()
Definition beammap.h:1162
Eigen::MatrixXd params
Definition beammap.h:43
auto run_timestream(KidsProc &)
Definition beammap.h:367
Definition engine.h:155
Eigen::Index hwpr_start_indices
Definition engine.h:208
void obsnum_setup()
Definition engine.h:353
std::string obsnum_dir_name
Definition engine.h:184
void write_stats()
Definition engine.h:2251
std::string parallel_policy
Definition engine.h:196
void write_map_summary(map_buffer_t &)
Definition engine.h:1584
std::map< std::string, Eigen::VectorXd > pointing_offsets_arcsec
Definition engine.h:241
std::vector< fitsIO< file_type_enum::write_fits, CCfits::ExtHDU * > > filtered_coadd_fits_io_vec
Definition engine.h:251
std::vector< Eigen::VectorXd > nw_times
Definition engine.h:166
std::vector< fitsIO< file_type_enum::write_fits, CCfits::ExtHDU * > > noise_fits_io_vec
Definition engine.h:246
std::vector< fitsIO< file_type_enum::write_fits, CCfits::ExtHDU * > > coadd_noise_fits_io_vec
Definition engine.h:250
std::string tod_output_type
Definition engine.h:223
std::string coadd_dir_name
Definition engine.h:184
std::vector< fitsIO< file_type_enum::write_fits, CCfits::ExtHDU * > > filtered_noise_fits_io_vec
Definition engine.h:247
bool verbose_mode
Definition engine.h:172
int n_maps
Definition engine.h:229
Eigen::VectorXI arrays_to_maps
Definition engine.h:232
std::vector< fitsIO< file_type_enum::write_fits, CCfits::ExtHDU * > > filtered_fits_io_vec
Definition engine.h:247
int n_scans_done
Definition engine.h:199
void add_tod_header()
Definition engine.h:1054
Eigen::VectorXd t_common
Definition engine.h:164
std::vector< fitsIO< file_type_enum::write_fits, CCfits::ExtHDU * > > filtered_coadd_noise_fits_io_vec
Definition engine.h:251
void add_phdu(fits_io_type &, map_buffer_t &, Eigen::Index)
Definition engine.h:1737
std::string redu_type
Definition engine.h:214
std::vector< Eigen::VectorXi > masks
Definition engine.h:165
std::string map_grouping
Definition engine.h:226
std::shared_ptr< spdlog::logger > logger
Definition engine.h:161
std::map< std::string, std::string > tod_filename
Definition engine.h:187
std::string tod_type
Definition engine.h:211
std::string obsnum
Definition engine.h:217
std::vector< fitsIO< file_type_enum::write_fits, CCfits::ExtHDU * > > coadd_fits_io_vec
Definition engine.h:250
std::vector< std::string > date_obs
Definition engine.h:169
void write_chunk_summary(TCData< tc_t, Eigen::MatrixXd > &)
Definition engine.h:1529
std::vector< fitsIO< file_type_enum::write_fits, CCfits::ExtHDU * > > fits_io_vec
Definition engine.h:246
int n_threads
Definition engine.h:193
std::string map_method
Definition engine.h:226
void write_maps(fits_io_type &, fits_io_type &, map_buffer_t &, Eigen::Index)
Definition engine.h:2030
Eigen::Index n_nws
Definition calib.h:50
Eigen::Index n_arrays
Definition calib.h:53
std::map< Eigen::Index, std::tuple< Eigen::Index, Eigen::Index > > array_limits
Definition calib.h:56
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
YAML::Node apt_meta
Definition calib.h:38
std::map< std::string, Eigen::VectorXd > apt
Definition calib.h:22
std::map< std::string, std::string > apt_header_description
Definition calib.h:151
Eigen::VectorXd hwpr_angle
Definition calib.h:24
std::map< Eigen::Index, std::tuple< Eigen::Index, Eigen::Index > > nw_limits
Definition calib.h:56
std::vector< std::string > apt_header_keys
Definition calib.h:81
Eigen::Index n_dets
Definition calib.h:47
void setup()
Definition calib.cpp:245
void calc_stats(timestream::TCData< tcdata_kind, Eigen::MatrixXd > &)
Definition diagnostics.h:61
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 source_name
Definition telescope.h:23
double d_fsmp
Definition telescope.h:42
double tau_225_GHz
Definition telescope.h:45
Eigen::MatrixXI scan_indices
Definition telescope.h:51
std::string project_id
Definition telescope.h:23
std::string pixel_axes
Definition telescope.h:59
std::map< std::string, Eigen::VectorXd > tel_header
Definition telescope.h:56
int n_params
Definition fitting.h:23
@ map
Definition toltec_io.h:22
std::map< Eigen::Index, std::string > array_name_map
Definition toltec_io.h:54
@ apt
Definition toltec_io.h:15
@ 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
void populate_maps_jinc(TCData< TCDataKind::PTC, Eigen::MatrixXd > &, map_buffer_t &, map_buffer_t &, Eigen::DenseBase< Derived > &, std::string &, apt_t &, double, bool, bool)
Definition jinc_mm.h:211
void populate_maps_jinc_parallel(TCData< TCDataKind::PTC, Eigen::MatrixXd > &, map_buffer_t &, map_buffer_t &, Eigen::DenseBase< Derived > &, std::string &, apt_t &, double, bool, bool)
Definition jinc_mm.h:485
Definition map.h:45
void calc_median_rms()
Definition map.cpp:416
void calc_map_hist()
Definition map.cpp:346
Eigen::Index n_rows
Definition map.h:65
std::string parallel_policy
Definition map.h:59
Eigen::Index n_cols
Definition map.h:65
void calc_map_psd()
Definition map.cpp:269
std::vector< Eigen::MatrixXd > coverage
Definition map.h:78
std::vector< Eigen::MatrixXd > signal
Definition map.h:78
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 normalize_maps()
Definition map.cpp:126
std::vector< Eigen::Tensor< double, 3 > > noise
Definition map.h:81
double pixel_size_rad
Definition map.h:69
void populate_maps_naive(TCData< TCDataKind::PTC, Eigen::MatrixXd > &, map_buffer_t &, map_buffer_t &, Eigen::DenseBase< Derived > &, std::string &, apt_t &, double, bool, bool)
Definition naive_mm.h:95
void populate_maps_naive_parallel(TCData< TCDataKind::PTC, Eigen::MatrixXd > &, map_buffer_t &, map_buffer_t &, Eigen::DenseBase< Derived > &, std::string &, apt_t &, double, bool, bool)
Definition naive_mm.h:439
auto calc_tau(Eigen::DenseBase< Derived > &, double)
Definition calibrate.h:136
Eigen::Index n_terms
Definition filter.h:19
Definition ptcproc.h:19
void append_to_netcdf(TCData< TCDataKind::PTC, Eigen::MatrixXd > &, std::string, std::string, std::string &, pointing_offset_t &, calib_t &)
Definition ptcproc.h:498
void calc_weights(TCData< TCDataKind::PTC, Eigen::MatrixXd > &, apt_type &, tel_type &)
Definition ptcproc.h:300
void run(TCData< TCDataKind::PTC, Eigen::MatrixXd > &, TCData< TCDataKind::PTC, Eigen::MatrixXd > &, calib_type &, std::string, std::string)
Definition ptcproc.h:173
auto reset_weights(TCData< TCDataKind::PTC, Eigen::MatrixXd > &, calib_t &, std::string)
Definition ptcproc.h:383
Definition rtcproc.h:18
void remove_flagged_dets(TCData< TCDataKind::PTC, Eigen::MatrixXd > &, apt_t &)
Definition rtcproc.h:493
bool run_polarization
Definition rtcproc.h:23
bool run_kernel
Definition rtcproc.h:24
timestream::Filter filter
Definition rtcproc.h:35
bool run_tod_filter
Definition rtcproc.h:26
auto remove_nearby_tones(TCData< TCDataKind::PTC, Eigen::MatrixXd > &, calib_t &, std::string)
Definition rtcproc.h:516
bool run_extinction
Definition rtcproc.h:29
timestream::Calibration calibration
Definition rtcproc.h:37
void append_to_netcdf(TCData< TCDataKind::PTC, Eigen::MatrixXd > &, std::string, std::string, std::string &, pointing_offset_t &, calib_t &)
Definition rtcproc.h:553
auto run(TCData< TCDataKind::RTC, Eigen::MatrixXd > &, TCData< TCDataKind::PTC, Eigen::MatrixXd > &, calib_t &, telescope_t &, double, std::string)
Definition rtcproc.h:261
void map_to_tod(mb_t &, TCData< tcdata_t, Eigen::MatrixXd > &, calib_t &, Eigen::DenseBase< Derived > &, std::string, std::string)
Definition timestream.h:624
auto remove_bad_dets(TCData< tcdata_t, Eigen::MatrixXd > &, calib_t &, std::string)
Definition timestream.h:703
@ Map
Definition timestream.h:271
@ Gaussian
Definition timestream.h:267
bool run_fruit_loops
Definition timestream.h:279
void add_gaussian(TCData< tcdata_t, Eigen::MatrixXd > &, Eigen::DenseBase< Derived > &, std::string &, std::string &, apt_t &, double, Eigen::Index, Eigen::Index)
Definition timestream.h:838
#define STD_TO_FWHM
Definition constants.h:45
#define mJY_ASEC_to_MJY_SR
Definition constants.h:54
#define RAD_TO_ASEC
Definition constants.h:36
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 double unix_to_modified_julian_date(double unix_time)
Definition utils.h:170
auto calc_det_pointing(tel_data_t &tel_data, double az_off, double el_off, const std::string pixel_axes, pointing_offset_t &pointing_offsets, const std::string map_grouping)
Definition pointing.h:11
auto tangent_to_abs(Eigen::DenseBase< Derived > &lat, Eigen::DenseBase< Derived > &lon, const double cra, const double cdec)
Definition pointing.h:88
static const std::string current_date_time()
Definition utils.h:91
auto calc_std_dev(Eigen::DenseBase< DerivedA > &data)
Definition utils.h:336
@ FilteredObs
Definition map.h:19
@ FilteredCoadd
Definition map.h:21
@ RawObs
Definition map.h:18
@ RawCoadd
Definition map.h:20
@ Good
Definition timestream.h:27
void calc_sensitivity(std::vector< tc_t > &ptcs, Eigen::DenseBase< DerivedA > &sensitivities, Eigen::DenseBase< DerivedB > &noisefluxes, double fsmp, Eigen::Index det, internal::Interval< double > freqrange={3., 5.})
Definition sensitivity.h:259
The raw obs struct This represents a single observation that contains a set of data items and calibra...
Definition io.h:50
bool beammap_subtract_reference
Definition engine.h:133
std::map< std::string, double > lower_sig2noise
Definition engine.h:148
std::map< std::string, double > beammap_fluxes_mJy_beam
Definition engine.h:123
Eigen::Index beammap_reference_det
Definition engine.h:136
std::map< std::string, double > upper_fwhm_arcsec
Definition engine.h:148
std::map< std::string, double > max_dist_arcsec
Definition engine.h:149
double lower_sens_factor
Definition engine.h:152
std::map< std::string, double > lower_fwhm_arcsec
Definition engine.h:148
int beammap_tod_output_iter
Definition engine.h:142
std::map< std::string, double > beammap_fluxes_MJy_Sr
Definition engine.h:124
bool beammap_derotate
Definition engine.h:139
std::map< std::string, double > upper_sig2noise
Definition engine.h:149
double upper_sens_factor
Definition engine.h:152
mapmaking::NaiveMapmaker naive_mm
Definition engine.h:109
engine::Telescope telescope
Definition engine.h:96
engine_utils::mapFitter map_fitter
Definition engine.h:99
engine::Diagnostics diagnostics
Definition engine.h:98
engine_utils::toltecIO toltec_io
Definition engine.h:97
timestream::PTCProc ptcproc
Definition engine.h:105
mapmaking::JincMapmaker jinc_mm
Definition engine.h:110
mapmaking::MapBuffer omb
Definition engine.h:108
mapmaking::MapBuffer cmb
Definition engine.h:108
timestream::RTCProc rtcproc
Definition engine.h:102
engine::Calib calib
Definition engine.h:95
bool run_noise
Definition engine.h:86
bool interp_over_gaps
Definition engine.h:73
bool run_tod_output
Definition engine.h:81
bool run_mapmaking
Definition engine.h:84
TC data class.
Definition timestream.h:55