Citlali
Loading...
Searching...
No Matches
todproc.h
Go to the documentation of this file.
1#pragma once
2
3#include <unsupported/Eigen/CXX11/Tensor>
4
5#include <tula/eigen.h>
6
8#include <tula/algorithm/mlinterp/mlinterp.hpp>
9
13
14namespace fs = std::filesystem;
15
17 template <typename OStream>
18 friend OStream &operator<<(OStream &os, const DummyEngine &e) {
19 return os << fmt::format("DummyEngine()");
20 }
21};
22
28template <class EngineType>
29struct TimeOrderedDataProc : ConfigMapper<TimeOrderedDataProc<EngineType>> {
31 using config_t = typename Base::config_t;
32 using Engine = EngineType;
33 using scanindicies_t = Eigen::MatrixXI;
34 using map_extent_t = std::vector<int>;
35 using map_coord_t = std::vector<Eigen::VectorXd>;
36 using map_count_t = std::size_t;
37 using array_indices_t = std::vector<std::tuple<Eigen::Index, Eigen::Index>>;
38 using det_indices_t = std::vector<std::tuple<Eigen::Index, Eigen::Index>>;
39
40 // get logger
41 std::shared_ptr<spdlog::logger> logger = spdlog::get("citlali_logger");
42
43 TimeOrderedDataProc(config_t config) : Base{std::move(config)} {}
44
45 // check if config file has nodes
46 static auto check_config(const config_t &config)
47 -> std::optional<std::string> {
48 // get logger
49 std::shared_ptr<spdlog::logger> logger = spdlog::get("citlali_logger");
50
51 std::vector<std::string> missing_keys;
52 logger->info("check TOD proc config\n{}", config);
53 // check for runtime config node
54 if (!config.has("runtime")) {
55 missing_keys.push_back("runtime");
56 }
57 // check for timestream config node
58 if (!config.has("timestream")) {
59 missing_keys.push_back("timestream");
60 }
61 // check for mapmaking config node
62 if (!config.has("mapmaking")) {
63 missing_keys.push_back("mapmaking");
64 }
65 // check for beammap config node
66 if (!config.has("beammap")) {
67 missing_keys.push_back("beammap");
68 }
69 // check for coadd config node
70 if (!config.has("coadd")) {
71 missing_keys.push_back("coadd");
72 }
73 // check for noise map config node
74 if (!config.has("noise_maps")) {
75 missing_keys.push_back("noise_maps");
76 }
77 // check for post processing config node
78 if (!config.has("post_processing")) {
79 missing_keys.push_back("post_processing");
80 }
81 if (missing_keys.empty()) {
82 return std::nullopt;
83 }
84 return fmt::format("invalid or missing keys={}", missing_keys);
85 }
86
87 // create output FITS files (does not populate)
89 // get apt from raw data files (beammapping)
90 void get_apt_from_files(const RawObs &rawobs);
91 // get tone frequencies from raw files
92 void get_tone_freqs_from_files(const RawObs &rawobs);
93 // get adc snap data from raw files
94 void get_adc_snap_from_files(const RawObs &rawobs);
95 // create output directories
96 void create_output_dir();
97 // count up detectors from input files and check for mismatch with apt
98 void check_inputs(const RawObs &rawobs);
99 // align networks and hwpr vectors in time
100 void align_timestreams(const RawObs &rawobs);
101 // updated alignment of networks and hwpr vectors in time
102 void align_timestreams_2(const RawObs &rawobs);
103 // interpolate pointing vectors
104 void interp_pointing();
105 // calculate number of maps
106 void calc_map_num();
107 // calculate size of omb maps
108 void calc_omb_size(std::vector<map_extent_t> &, std::vector<map_coord_t> &);
109 // allocate observation maps
111 // calculate size of cmb maps
112 void calc_cmb_size(std::vector<map_coord_t> &);
113 // allocate coadded maps
114 void allocate_cmb();
115 // allocate noise maps
116 template<class map_buffer_t>
117 void allocate_nmb(map_buffer_t &);
118 // coadd omb into cmb
119 void coadd();
120 // make index files
121 void make_index_file(std::string);
122
123 // TODO fix the const correctness
124 Engine &engine() { return m_engine; }
125
126 const Engine &engine() const { return m_engine; }
127
128 template <typename OStream>
129 friend OStream &operator<<(OStream &os,
130 const TimeOrderedDataProc &todproc) {
131 return os << fmt::format("TimeOrderedDataProc(engine={})",
132 todproc.engine());
133 }
134
135private:
137};
138
139// make apt table from raw files instead of an ecsv table
140template <class EngineType>
142 using namespace netCDF;
143 using namespace netCDF::exceptions;
144
145 // nw names
146 std::vector<Eigen::Index> interfaces;
147
148 // total number of detectors
149 Eigen::Index n_dets = 0;
150 // detector, nw and array names for each network
151 std::vector<Eigen::Index> dets, nws, arrays;
152 // loop through input files
153 for (const RawObs::DataItem &data_item : rawobs.kidsdata()) {
154 try {
155 // load data file
156 NcFile fo(data_item.filepath(), NcFile::read);
157 auto vars = fo.getVars();
158
159 // get the interface
160 interfaces.push_back(std::stoi(data_item.interface().substr(6)));
161
162 // add the current file's number of dets to the total
163 n_dets += vars.find("Data.Toltec.Is")->second.getDim(1).getSize();
164
165 // get the number of dets in file
166 dets.push_back(vars.find("Data.Toltec.Is")->second.getDim(1).getSize());
167 // get the nw from interface
168 nws.push_back(interfaces.back());
169 // get the array from the interface
170 arrays.push_back(engine().toltec_io.nw_to_array_map[interfaces.back()]);
171
172 fo.close();
173
174 } catch (NcException &e) {
175 logger->error("{}", e.what());
176 throw DataIOError{fmt::format(
177 "failed to load data from netCDF file {}", data_item.filepath())};
178 }
179 }
180
181 // explicitly clear the apt
182 engine().calib.apt.clear();
183
184 // resize the apt vectors
185 for (auto const& key : engine().calib.apt_header_keys) {
186 if (key=="x_t" || key=="y_t") {
187 engine().calib.apt[key].setZero(n_dets);
188 }
189 else {
190 engine().calib.apt[key].setOnes(n_dets);
191 }
192 }
193
194 // set all flags to good
195 engine().calib.apt["flag"].setZero(n_dets);
196
197 // add the nws and arrays to the apt table
198 Eigen::Index j = 0;
199 for (Eigen::Index i=0; i<nws.size(); ++i) {
200 engine().calib.apt["nw"].segment(j,dets[i]).setConstant(nws[i]);
201 engine().calib.apt["array"].segment(j,dets[i]).setConstant(arrays[i]);
202
203 j = j + dets[i];
204 }
205
206 // set uids
207 engine().calib.apt["uid"] = Eigen::VectorXd::LinSpaced(n_dets,0,n_dets-1);
208
209 // setup nws, arrays, etc.
210 engine().calib.setup();
211
212 // filepath
213 engine().calib.apt_filepath = "internally generated for beammap";
214}
215
216template <class EngineType>
218 using namespace netCDF;
219 using namespace netCDF::exceptions;
220
221 // tone frquencies for each network
222 std::map<Eigen::Index,Eigen::MatrixXd> tone_freqs;
223
224 // nw names
225 std::vector<Eigen::Index> interfaces;
226
227 // total number of detectors
228 Eigen::Index n_dets = 0;
229 // detector, nw and array names for each network
230 std::vector<Eigen::Index> dets, nws, arrays;
231 // loop through input files
232 for (const RawObs::DataItem &data_item : rawobs.kidsdata()) {
233 try {
234 // load data file
235 NcFile fo(data_item.filepath(), NcFile::read);
236 auto vars = fo.getVars();
237
238 // get the interface
239 interfaces.push_back(std::stoi(data_item.interface().substr(6)));
240
241 // add the current file's number of dets to the total
242 n_dets += vars.find("Data.Toltec.Is")->second.getDim(1).getSize();
243
244 // dimension of tone freqs is (n_sweeps, n_tones)
245 Eigen::Index n_sweeps = vars.find("Header.Toltec.ToneFreq")->second.getDim(0).getSize();
246
247 // get local oscillator frequency
248 double lo_freq;
249 vars.find("Header.Toltec.LoCenterFreq")->second.getVar(&lo_freq);
250
251 // get tone_freqs for interface
252 tone_freqs[interfaces.back()].resize(vars.find("Header.Toltec.ToneFreq")->second.getDim(1).getSize(),n_sweeps);
253 vars.find("Header.Toltec.ToneFreq")->second.getVar(tone_freqs[interfaces.back()].data());
254
255 // add local oscillator freq
256 tone_freqs[interfaces.back()] = tone_freqs[interfaces.back()].array() + lo_freq;
257
258 fo.close();
259
260 } catch (NcException &e) {
261 logger->error("{}", e.what());
262 throw DataIOError{fmt::format(
263 "failed to load data from netCDF file {}", data_item.filepath())};
264 }
265 }
266
267 engine().calib.apt["tone_freq"].resize(engine().calib.n_dets);
268
269 // add the nws and arrays to the apt table
270 Eigen::Index j = 0;
271 for (Eigen::Index i=0; i<engine().calib.nws.size(); ++i) {
272 engine().calib.apt["tone_freq"].segment(j,tone_freqs[interfaces[i]].size()) = tone_freqs[interfaces[i]];
273
274 j = j + tone_freqs[interfaces[i]].size();
275 }
276
277 if (!engine().telescope.sim_obs) {
278 /* find duplicates */
279
280 // frequency separation
281 Eigen::VectorXd dfreq(engine().calib.n_dets);
282 dfreq(0) = engine().calib.apt["tone_freq"](1) - engine().calib.apt["tone_freq"](0);
283
284 // loop through tone freqs and find distance
285 for (Eigen::Index i=1; i<engine().calib.apt["tone_freq"].size()-1; ++i) {
286 dfreq(i) = std::min(abs(engine().calib.apt["tone_freq"](i) - engine().calib.apt["tone_freq"](i-1)),
287 abs(engine().calib.apt["tone_freq"](i+1) - engine().calib.apt["tone_freq"](i)));
288 }
289 // get last distance
290 dfreq(dfreq.size()-1) = abs(engine().calib.apt["tone_freq"](dfreq.size()-1)-engine().calib.apt["tone_freq"](dfreq.size()-2));
291
292 // number of nearby tones found
293 int n_nearby_tones = 0;
294
295 // store duplicates
296 engine().calib.apt["duplicate_tone"].setZero(engine().calib.n_dets);
297
298 // loop through flag columns
299 for (Eigen::Index i=0; i<engine().calib.n_dets; ++i) {
300 // if closer than freq separation limit and unflagged, flag it
301 if (dfreq(i) < engine().rtcproc.delta_f_min_Hz) {
302 engine().calib.apt["duplicate_tone"](i) = 1;
303 n_nearby_tones++;
304 }
305 }
306 logger->info("{} nearby tones found. these will be flagged.",n_nearby_tones);
307 }
308}
309
310template <class EngineType>
312 using namespace netCDF;
313 using namespace netCDF::exceptions;
314
315 // explicitly clear adc vector
316 engine().diagnostics.adc_snap_data.clear();
317
318 // loop through input files
319 for (const RawObs::DataItem &data_item : rawobs.kidsdata()) {
320 try {
321 // load data file
322 NcFile fo(data_item.filepath(), NcFile::read);
323 auto vars = fo.getVars();
324
325 // dimension 0 of adc data
326 Eigen::Index adcSnapDim = vars.find("Header.Toltec.AdcSnapData")->second.getDim(0).getSize();
327 // dimension 1 of adc data
328 Eigen::Index adcSnapDataDim = vars.find("Header.Toltec.AdcSnapData")->second.getDim(1).getSize();
329
330 // matrix to hold adc data for current file
331 Eigen::Matrix<short,Eigen::Dynamic, Eigen::Dynamic> adcsnap(adcSnapDataDim,adcSnapDim);
332 // load adc data
333 vars.find("Header.Toltec.AdcSnapData")->second.getVar(adcsnap.data());
334 // append to vector of adc data
335 engine().diagnostics.adc_snap_data.push_back(adcsnap);
336
337 fo.close();
338
339 } catch (NcException &e) {
340 logger->warn("{} adc data not found",data_item.filepath());
341 }
342 }
343}
344
345// create output directories
346template <class EngineType>
348 // redu subdir
349 engine().redu_dir_name = "";
350
351 // create reduction subdir
352 if (engine().use_subdir) {
353 // redu number
354 engine().redu_dir_num = 0;
355
356 std::stringstream ss_redu_dir_num;
357 // add leading zero to redu_dir_num (i.e., '00', '01',...)
358 ss_redu_dir_num << std::setfill('0') << std::setw(2) << engine().redu_dir_num;
359
360 // create redu dir name ('redu00', 'redu01',...)
361 std::string redu_dir_name = "redu" + ss_redu_dir_num.str();
362
363 // iteratively check if current subdir with current redu number exists
364 while (fs::exists(fs::status(engine().output_dir + "/" + redu_dir_name))) {
365 // increment redu number if subdir exists
366 engine().redu_dir_num++;
367 std::stringstream ss_redu_dir_num_i;
368 ss_redu_dir_num_i << std::setfill('0') << std::setw(2) << engine().redu_dir_num;
369 redu_dir_name = "redu" + ss_redu_dir_num_i.str();
370 }
371
372 // final redu dir name is output directory from config + /reduNN
373 engine().redu_dir_name = engine().output_dir + "/" + redu_dir_name;
374
375 // create redu dir directory
376 fs::create_directories(engine().redu_dir_name);
377 }
378 else {
379 engine().redu_dir_name = engine().output_dir + "/";
380 }
381
382 // coadded subdir
383 if (engine().run_coadd) {
384 engine().coadd_dir_name = engine().redu_dir_name + "/coadded/";
385 // coadded raw subdir
386 if (!fs::exists(fs::status(engine().coadd_dir_name + "raw/"))) {
387 fs::create_directories(engine().coadd_dir_name + "raw/");
388 }
389 else {
390 logger->warn("directory {} already exists", engine().coadd_dir_name + "raw/");
391 }
392 // if map filtering is requested
393 if (engine().run_map_filter) {
394 // coadded filtered subdir
395 if (!fs::exists(fs::status(engine().coadd_dir_name + "filtered/"))) {
396 fs::create_directories(engine().coadd_dir_name + "filtered/");
397 }
398 else {
399 logger->warn("directory {} already exists", engine().coadd_dir_name + "filtered/");
400 }
401 }
402 }
403}
404
405template <class EngineType>
407 using namespace netCDF;
408 using namespace netCDF::exceptions;
409
410 Eigen::Index n_dets = 0;
411
412 // loop through input files
413 for (const RawObs::DataItem &data_item : rawobs.kidsdata()) {
414 try {
415 // load data file
416 NcFile fo(data_item.filepath(), NcFile::read);
417 auto vars = fo.getVars();
418 // get number of dets from data and add to global value
419 n_dets += vars.find("Data.Toltec.Is")->second.getDim(1).getSize();
420
421 fo.close();
422
423 } catch (NcException &e) {
424 logger->error("{}", e.what());
425 throw DataIOError{fmt::format(
426 "failed to load data from netCDF file {}", data_item.filepath())};
427 }
428 }
429
430 // check if number of detectors in apt file is equal to those in files
431 if (n_dets != engine().calib.n_dets) {
432 logger->error("number of detectors in data files and apt file do not match");
433 std::exit(EXIT_FAILURE);
434 }
435}
436
437// align tod with telescope
438template <class EngineType>
440 using namespace netCDF;
441 using namespace netCDF::exceptions;
442
443 // clear start and end indices for each observation
444 engine().start_indices.clear();
445 engine().end_indices.clear();
446
447 // clear gaps
448 engine().gaps.clear();
449
450 // vector of network times
451 std::vector<Eigen::VectorXd> nw_ts;
452 // start and end times
453 std::vector<double> nw_t0, nw_tn;
454
455 // maximum start time
456 double max_t0 = -99;
457
458 // minimum end time
459 double min_tn = std::numeric_limits<double>::max();
460 // indices of max start time and min end time
461 Eigen::Index max_t0_i, min_tn_i;
462
463 // set network
464 Eigen::Index nw = 0;
465 // sample rate
466 double fsmp = -1;
467
468 // loop through input files
469 for (const RawObs::DataItem &data_item : rawobs.kidsdata()) {
470 try {
471 // load data file
472 NcFile fo(data_item.filepath(), NcFile::read);
473 auto vars = fo.getVars();
474
475 // get roach index for offsets
476 int roach_index;
477 vars.find("Header.Toltec.RoachIndex")->second.getVar(&roach_index);
478
479 // get sample rate
480 double fsmp_roach;
481 vars.find("Header.Toltec.SampleFreq")->second.getVar(&fsmp_roach);
482
483 // check if sample rate is the same and exit if not
484 if (fsmp!=-1 && fsmp_roach!=fsmp) {
485 logger->error("mismatched sample rate in toltec{}",roach_index);
486 std::exit(EXIT_FAILURE);
487 }
488 else {
489 fsmp = fsmp_roach;
490 }
491
492 // get dimensions for time matrix
493 Eigen::Index n_pts = vars.find("Data.Toltec.Ts")->second.getDim(0).getSize();
494 Eigen::Index n_time = vars.find("Data.Toltec.Ts")->second.getDim(1).getSize();
495
496 // get time matrix
497 Eigen::MatrixXi ts(n_time,n_pts);
498 vars.find("Data.Toltec.Ts")->second.getVar(ts.data());
499
500 // transpose due to row-major order
501 ts.transposeInPlace();
502
503 // find gaps
504 int gaps = ((ts.block(1,3,n_pts,1).array() - ts.block(0,3,n_pts-1,1).array()).array() > 1).count();
505
506 // add gaps to engine map
507 if (gaps>0) {
508 engine().gaps["Toltec" + std::to_string(roach_index)] = gaps;
509 }
510
511 // get fpga frequency
512 double fpga_freq;
513 vars.find("Header.Toltec.FpgaFreq")->second.getVar(&fpga_freq);
514
515 // ClockTime (sec)
516 auto sec0 = ts.cast <double> ().col(0);
517 // ClockTimeNanoSec (nsec)
518 auto nsec0 = ts.cast <double> ().col(5);
519 // PpsCount (pps ticks)
520 auto pps = ts.cast <double> ().col(1);
521 // ClockCount (clock ticks)
522 auto msec = ts.cast <double> ().col(2)/fpga_freq;
523 // PacketCount (packet ticks)
524 auto count = ts.cast <double> ().col(3);
525 // PpsTime (clock ticks)
526 auto pps_msec = ts.cast <double> ().col(4)/fpga_freq;
527 // get start time
528 auto t0 = sec0 + nsec0*1e-9;
529
530 // shift start time (offset determined empirically)
531 int start_t = int(t0[0] - 0.5);
532 //int start_t = int(t0[0]);
533
534 // convert start time to double
535 double start_t_dbl = start_t;
536 // clock count - clock ticks
537 Eigen::VectorXd dt = msec - pps_msec;
538 // remove overflow due to int32
539 dt = (dt.array() < 0).select(msec.array() - pps_msec.array() + (pow(2.0,32)-1)/fpga_freq,msec - pps_msec);
540 // get network time and add offsets
541 nw_ts.push_back(start_t_dbl + pps.array() + dt.array() +
542 engine().interface_sync_offset["toltec"+std::to_string(roach_index)]);
543
544 // push back start time
545 nw_t0.push_back(nw_ts.back()[0]);
546
547 // push back end time
548 nw_tn.push_back(nw_ts.back()[n_pts - 1]);
549
550 // get global max start time and index
551 if (nw_t0.back() > max_t0) {
552 max_t0 = nw_t0.back();
553 max_t0_i = nw;
554 }
555
556 // get global min end time and index
557 if (nw_tn.back() < min_tn) {
558 min_tn = nw_tn.back();
559 min_tn_i = nw;
560 }
561
562 // increment nw
563 nw++;
564
565 fo.close();
566
567 } catch (NcException &e) {
568 logger->error("{}", e.what());
569 throw DataIOError{fmt::format(
570 "failed to load data from netCDF file {}", data_item.filepath())};
571 }
572 }
573
574 // get hwpr timing
575 if (engine().calib.run_hwpr) {
576 /*auto sec0 = engine().calib.hwpr_ts.template cast <double> ().col(0);
577 // ClockTimeNanoSec (nsec)
578 auto nsec0 = engine().calib.hwpr_ts.template cast <double> ().col(5);
579 // PpsCount (pps ticks)
580 auto pps = engine().calib.hwpr_ts.template cast <double> ().col(1);
581 // ClockCount (clock ticks)
582 auto msec = engine().calib.hwpr_ts.template cast <double> ().col(2)/engine().calib.hwpr_fpga_freq;
583 // PacketCount (packet ticks)
584 auto count = engine().calib.hwpr_ts.template cast <double> ().col(3);
585 // PpsTime (clock ticks)
586 auto pps_msec = engine().calib.hwpr_ts.template cast <double> ().col(4)/engine().calib.hwpr_fpga_freq;
587 // get start time
588 auto t0 = sec0 + nsec0*1e-9;
589
590 // shift start time
591 int start_t = int(t0[0] - 0.5);
592 //int start_t = int(t0[0]);
593
594 // convert start time to double
595 double start_t_dbl = start_t;
596
597 Eigen::VectorXd dt = msec - pps_msec;
598
599 // remove overflow due to int32
600 dt = (dt.array() < 0).select(msec.array() - pps_msec.array() + (pow(2.0,32)-1)/engine().calib.hwpr_,msec - pps_msec);
601
602 // get network time and add offsets
603 engine().calib.hwpr_recvt = start_t_dbl + pps.array() + dt.array() + engine().interface_sync_offset["hwpr"];
604 */
605
606 // if hwpr init time is larger than max start time, replace global max start time
607 Eigen::Index hwpr_ts_n_pts = engine().calib.hwpr_recvt.size();
608 if (engine().calib.hwpr_recvt(0) > max_t0) {
609 max_t0 = engine().calib.hwpr_recvt(0);
610 }
611
612 // if hwpr init time is smaller than min end time, replace global min end time
613 if (engine().calib.hwpr_recvt(hwpr_ts_n_pts - 1) < min_tn) {
614 min_tn = engine().calib.hwpr_recvt(hwpr_ts_n_pts - 1);
615 }
616 }
617
618 // size of smallest data time vector
619 Eigen::Index min_size = nw_ts[0].size();
620
621 // loop through time vectors and get the smallest
622 for (Eigen::Index i=0; i<nw_ts.size(); ++i) {
623 // find start index that is larger than max start
624 Eigen::Index si, ei;
625 // find index closest to max start time
626 auto s = (abs(nw_ts[i].array() - max_t0)).minCoeff(&si);
627
628 // if closest index is smaller than max start time
629 // incrememnt index until it is larger or equal
630 while (nw_ts[i][si] < max_t0) {
631 si++;
632 }
633 // pushback start index on start index vector
634 engine().start_indices.push_back(si);
635
636 // find end index that is smaller than min end
637 auto e = (abs(nw_ts[i].array() - min_tn)).minCoeff(&ei);
638 // if closest index is larger than min end time
639 // incrememnt index until it is smaller or equal
640 while (nw_ts[i][ei] > min_tn) {
641 ei--;
642 }
643 // pushback end index on end index vector
644 engine().end_indices.push_back(ei);
645 }
646
647 // get min size
648 for (Eigen::Index i=0; i<nw_ts.size(); ++i) {
649 // start indices
650 auto si = engine().start_indices[i];
651 // end indices
652 auto ei = engine().end_indices[i];
653
654 // if smallest length, update min_size
655 if ((ei - si + 1) < min_size) {
656 min_size = ei - si + 1;
657 }
658 }
659
660 // if hwpr requested
661 if (engine().calib.run_hwpr) {
662 Eigen::Index si, ei;
663 // find start index that is larger than max start for hwpr
664 auto s = (abs(engine().calib.hwpr_recvt.array() - max_t0)).minCoeff(&si);
665
666 // if closest index is smaller than max start time
667 // incrememnt index until it is larger or equal
668 while (engine().calib.hwpr_recvt(si) < max_t0) {
669 si++;
670 }
671 // pushback start index on hwpr start index vector
672 engine().hwpr_start_indices = si;
673
674 // find end index that is smaller than min end for hwpr
675 auto e = (abs(engine().calib.hwpr_recvt.array() - min_tn)).minCoeff(&ei);
676 // if closest index is larger than min end time
677 // incrememnt index until it is smaller or equal
678 while (engine().calib.hwpr_recvt(ei) > min_tn) {
679 ei--;
680 }
681 // pushback end index on hwpr end index vector
682 engine().hwpr_end_indices = ei;
683
684 // update min_size for all time vectors if hwpr data is shorter (data and hwpr)
685 if ((ei - si + 1) < min_size) {
686 min_size = ei - si + 1;
687 }
688 }
689
690 // size of telescope data
691 Eigen::Matrix<Eigen::Index,1,1> nd;
692 nd << engine().telescope.tel_data["TelTime"].size();
693
694 // shortest data time vector
695 Eigen::VectorXd xi = nw_ts[max_t0_i].head(min_size);
696
697 // interpolate telescope data
698 for (const auto &tel_it : engine().telescope.tel_data) {
699 if (tel_it.first !="TelTime") {
700 // telescope vector to interpolate
701 Eigen::VectorXd yd = engine().telescope.tel_data[tel_it.first];
702 // vector to store interpolated outputs in
703 Eigen::VectorXd yi(min_size);
704
705 mlinterp::interp(nd.data(), min_size, // nd, ni
706 yd.data(), yi.data(), // yd, yi
707 engine().telescope.tel_data["TelTime"].data(), xi.data()); // xd, xi
708
709 // move back into tel_data vector
710 engine().telescope.tel_data[tel_it.first] = std::move(yi);
711 }
712 }
713
714 // replace telescope time vectors
715 engine().telescope.tel_data["TelTime"] = xi;
716 engine().telescope.tel_data["TelUTC"] = xi;
717
718 // interpolate hwpr data
719 if (engine().calib.run_hwpr) {
720 Eigen::VectorXd yd = engine().calib.hwpr_angle;
721 // vector to store interpolated outputs in
722 Eigen::VectorXd yi(min_size);
723 mlinterp::interp(nd.data(), min_size, // nd, ni
724 yd.data(), yi.data(), // yd, yi
725 engine().calib.hwpr_recvt.data(), xi.data()); // xd, xi
726
727 // move back into hwpr angle
728 engine().calib.hwpr_angle = std::move(yi);
729 }
730}
731
732// upgraded alignment of tod with telescope
733template <class EngineType>
735
736template <class EngineType>
738 // how many offsets in config file
739 Eigen::Index n_offsets = engine().pointing_offsets_arcsec["az"].size();
740
741 // keys for pointing offsets
742 std::vector<std::string> altaz_keys = {"alt", "az"};
743
744 for (const auto &key: altaz_keys) {
745 // if only one value given
746 if (n_offsets==1) {
747 double offset = engine().pointing_offsets_arcsec[key](0);
748 engine().pointing_offsets_arcsec[key].resize(engine().telescope.tel_data["TelTime"].size());
749 engine().pointing_offsets_arcsec[key].setConstant(offset);
750 }
751 else if (n_offsets==2) {
752 // size of telescope data
753 Eigen::Index ni = engine().telescope.tel_data["TelTime"].size();
754
755 // size of telescope data
756 Eigen::Matrix<Eigen::Index,1,1> nd;
757 nd << n_offsets;
758
759 // vector to store interpolation
760 Eigen::VectorXd yi(ni);
761
762 // start and end times of observation
763 Eigen::VectorXd xd(n_offsets);
764 // use start and end of current obs if julian dates not specified
765 if ((engine().pointing_offsets_modified_julian_date<=0).any()) {
766 xd << engine().telescope.tel_data["TelTime"](0), engine().telescope.tel_data["TelTime"](ni-1);
767 }
768 // else use specified modified julian dates, convert to julian dates, and calc unix time
769 else {
770 xd << engine_utils::modified_julian_date_to_unix(engine().pointing_offsets_modified_julian_date(0)),
771 engine_utils::modified_julian_date_to_unix(engine().pointing_offsets_modified_julian_date(1));
772
773 // make sure offsets are before and after the observation
774 if (xd(0) > engine().telescope.tel_data["TelTime"](0) || xd(1) < engine().telescope.tel_data["TelTime"](ni-1)) {
775 logger->error("MJD range is invalid");
776 std::exit(EXIT_FAILURE);
777 }
778 }
779
780 // interpolate offset onto time vector
781 mlinterp::interp(nd.data(), ni, // nd, ni
782 engine().pointing_offsets_arcsec[key].data(), yi.data(), // yd, yi
783 xd.data(), engine().telescope.tel_data["TelTime"].data()); // xd, xi
784
785 // overwrite pointing offsets
786 engine().pointing_offsets_arcsec[key] = yi;
787 }
788 else {
789 logger->error("only one or two values for altaz offsets are supported");
790 std::exit(EXIT_FAILURE);
791 }
792 }
793}
794
795// get map number
796template <class EngineType>
798 // auto map grouping
799 if (engine().map_grouping=="auto") {
800 // array map grouping for science and pointing
801 if ((engine().redu_type == "science") || (engine().redu_type == "pointing")) {
802 engine().map_grouping = "array";
803 }
804
805 // detector map grouping for beammaps
806 else if ((engine().redu_type == "beammap")) {
807 engine().map_grouping = "detector";
808 }
809 }
810
811 // overwrite map number for detectors
812 if (engine().map_grouping == "detector") {
813 engine().n_maps = engine().calib.n_dets;
814 }
815
816 // overwrite map number for networks
817 else if (engine().map_grouping == "nw") {
818 engine().n_maps = engine().calib.n_nws;
819 }
820
821 // overwrite map number for arrays
822 else if (engine().map_grouping == "array") {
823 engine().n_maps = engine().calib.n_arrays;
824 }
825
826 // overwrite map number for fg grouping
827 else if (engine().map_grouping == "fg") {
828 // there are potentially 4 fg's per array, so total number of maps is max 4 x n_arrays
829 engine().n_maps = engine().calib.fg.size()*engine().calib.n_arrays;
830 }
831
832 if (engine().rtcproc.run_polarization) {
833 // multiply by number of polarizations (stokes I + Q + U = 3)
834 engine().n_maps = engine().n_maps*engine().rtcproc.polarization.stokes_params.size();
835 }
836
837 // mapping from index in map vector to detector array index
838 // if stokes I array grouping with all arrays, this will be [0,1,2]
839 // if missing array 0, this will be [1,2]
840 engine().maps_to_arrays.resize(engine().n_maps);
841
842 // mapping from index in map vector to stokes parameter index (I=0, Q=1, U=2)
843 // if array grouping with all arrays this will be [0,0,0,1,1,2,2,2]
844 // and maps_to_arrays will be [0,1,2,0,1,2,0,1,2]
845 engine().maps_to_stokes.resize(engine().n_maps);
846
847 // mapping from array index to index in map vectors (reverse of maps_to_arrays)
848 // if stokes I array grouping with all arrays, this will also be [0,1,2]
849 // if missing array 0, this will be [0,1]
850 engine().arrays_to_maps.resize(engine().n_maps);
851
852 // array to hold mapping from group to detector array index
853 Eigen::VectorXI array_indices;
854
855 // detector gropuing
856 if (engine().map_grouping == "detector") {
857 // only do stokes I as Q and U don't make sense for detector grouping
858 // this is just a copy of the array indices from the apt
859 array_indices = engine().calib.apt["array"].template cast<Eigen::Index> ();
860 }
861
862 // array grouping
863 else if (engine().map_grouping == "array") {
864 // if all arrays are included this will be [0,1,2]
865 array_indices = engine().calib.arrays;
866 }
867
868 // network grouping
869 else if (engine().map_grouping == "nw") {
870 // if all nws/arrays are included this will be:
871 // [0,0,0,0,0,0,0,0,1,1,1,1,2,2]
872 // nws are ordered automatically when files are read in
873 array_indices.resize(engine().calib.nws.size());
874
875 // find all map from nw to arrays
876 for (Eigen::Index i=0; i<engine().calib.nws.size(); ++i) {
877 // get array for current nw
878 array_indices(i) = engine().toltec_io.nw_to_array_map[engine().calib.nws(i)];
879 }
880 }
881
882 // frequency grouping
883 else if (engine().map_grouping == "fg") {
884 // size of array indices is number of fg's x number of arrays
885 // if all fgs are included, this will be:
886 // [0,0,0,0,0,1,1,1,1,1,2,2,2,2,2]
887 // the order of the fgs will vary depending on the apt, but this is irrelevant
888 array_indices.resize(engine().calib.fg.size()*engine().calib.n_arrays);
889
890 // map from fg to array index
891 Eigen::Index j = 0;
892 // loop through arrays
893 for (Eigen::Index i=0; i<engine().calib.n_arrays; ++i) {
894 // append current array index to all elements within a segment of fg size
895 array_indices.segment(j,engine().calib.fg.size()).setConstant(engine().calib.arrays(i));
896 // increment by fg size
897 j = j + engine().calib.fg.size();
898 }
899 }
900
901 // copy array_indices into maps_to_arrays and maps_to_stokes for each stokes param
902 Eigen::Index j = 0;
903 // loop through stokes params
904 for (const auto &[stokes_index,stokes_param]: engine().rtcproc.polarization.stokes_params) {
905 // for each stokes param append all array indices in order
906 engine().maps_to_arrays.segment(j,array_indices.size()) = array_indices;
907 // for each stokes param append current stokes index
908 engine().maps_to_stokes.segment(j,array_indices.size()).setConstant(stokes_index);
909 // increment by array index size
910 j = j + array_indices.size();
911 }
912
913 // calculate detector array index to map index
914 Eigen::Index index = 0;
915 // start at map index 0
916 engine().arrays_to_maps(0) = index;
917 for (Eigen::Index i=1; i<engine().n_maps; ++i) {
918 // we move to the next map index when the array increments
919 if (engine().maps_to_arrays(i) > engine().maps_to_arrays(i-1)) {
920 index++;
921 }
922 // reset to first map index when we return the an earlier array
923 else if (engine().maps_to_arrays(i) < engine().maps_to_arrays(i-1)) {
924 index = 0;
925 }
926 engine().arrays_to_maps(i) = index;
927 }
928}
929
930// calculate map dimensions
931template <class EngineType>
932void TimeOrderedDataProc<EngineType>::calc_omb_size(std::vector<map_extent_t> &map_extents, std::vector<map_coord_t> &map_coords) {
933
934 // reference to map buffer
935 auto& omb = engine().omb;
936
937 // only run if manual map sizes have not been input
938 if ((engine().omb.wcs.naxis[0] <= 0) || (engine().omb.wcs.naxis[1] <= 0)) {
939 // matrix to store size limits
940 Eigen::MatrixXd det_lat_limits(engine().calib.n_dets, 2);
941 Eigen::MatrixXd det_lon_limits(engine().calib.n_dets, 2);
942 det_lat_limits.setZero();
943 det_lon_limits.setZero();
944
945 // placeholder vectors for grppi maps
946 std::vector<int> det_in_vec, det_out_vec;
947
948 // placeholder vectors for grppi loop
949 det_in_vec.resize(engine().calib.n_dets);
950 std::iota(det_in_vec.begin(), det_in_vec.end(), 0);
951 det_out_vec.resize(engine().calib.n_dets);
952
953 // get telescope meta data for current scan
954 std::map<std::string, Eigen::VectorXd> tel_data;
955
956 // pointing offsets
957 std::map<std::string, Eigen::VectorXd> pointing_offsets_arcsec;
958
959 // loop through scans
960 for (Eigen::Index i=0; i<engine().telescope.scan_indices.cols(); ++i) {
961 // lower scan index
962 auto si = engine().telescope.scan_indices(0,i);
963 // upper scan index
964 auto sl = engine().telescope.scan_indices(1,i) - engine().telescope.scan_indices(0,i) + 1;
965
966 for (auto const& x: engine().telescope.tel_data) {
967 tel_data[x.first] = engine().telescope.tel_data[x.first].segment(si,sl);
968 }
969
970 // get pointing offsets for current scan
971 pointing_offsets_arcsec["az"] = engine().pointing_offsets_arcsec["az"].segment(si,sl);
972 pointing_offsets_arcsec["alt"] = engine().pointing_offsets_arcsec["alt"].segment(si,sl);
973
974 // don't need to find the offsets if in detector mode
975 if (engine().map_grouping!="detector") {
976 // loop through detectors
977 grppi::map(tula::grppi_utils::dyn_ex(engine().parallel_policy), det_in_vec, det_out_vec, [&](auto j) {
978
979 // get pointing
980 auto [lat, lon] = engine_utils::calc_det_pointing(tel_data, engine().calib.apt["x_t"](j), engine().calib.apt["y_t"](j),
981 engine().telescope.pixel_axes, pointing_offsets_arcsec, engine().map_grouping);
982 // check for min and max
983 if (engine().calib.apt["flag"](j)==0) {
984 if (lat.minCoeff() < det_lat_limits(j,0)) {
985 det_lat_limits(j,0) = lat.minCoeff();
986 }
987 if (lat.maxCoeff() > det_lat_limits(j,1)) {
988 det_lat_limits(j,1) = lat.maxCoeff();
989 }
990 if (lon.minCoeff() < det_lon_limits(j,0)) {
991 det_lon_limits(j,0) = lon.minCoeff();
992 }
993 if (lon.maxCoeff() > det_lon_limits(j,1)) {
994 det_lon_limits(j,1) = lon.maxCoeff();
995 }
996 }
997 return 0;
998 });
999 }
1000 else {
1001 // calculate detector pointing for first detector only since offsets are zero
1002 auto [lat, lon] = engine_utils::calc_det_pointing(tel_data, 0., 0., engine().telescope.pixel_axes,
1003 pointing_offsets_arcsec, engine().map_grouping);
1004 if (lat.minCoeff() < det_lat_limits(0,0)) {
1005 det_lat_limits.col(0).setConstant(lat.minCoeff());
1006 }
1007 if (lat.maxCoeff() > det_lat_limits(0,1)) {
1008 det_lat_limits.col(1).setConstant(lat.maxCoeff());
1009 }
1010 if (lon.minCoeff() < det_lon_limits(0,0)) {
1011 det_lon_limits.col(0).setConstant(lon.minCoeff());
1012 }
1013 if (lon.maxCoeff() > det_lon_limits(0,1)) {
1014 det_lon_limits.col(1).setConstant(lon.maxCoeff());
1015 }
1016 }
1017 }
1018
1019 // get the global min and max
1020 double min_lat = det_lat_limits.col(0).minCoeff();
1021 double max_lat = det_lat_limits.col(1).maxCoeff();
1022 double min_lon = det_lon_limits.col(0).minCoeff();
1023 double max_lon = det_lon_limits.col(1).maxCoeff();
1024
1025 // calculate dimensions
1026 auto calc_map_dims = [&](double min_dim, double max_dim) {
1027 int min_pix = static_cast<int>(ceil(abs(min_dim / omb.pixel_size_rad)));
1028 int max_pix = static_cast<int>(ceil(abs(max_dim / omb.pixel_size_rad)));
1029 return 2 * std::max(min_pix, max_pix) + 1;
1030 };
1031
1032 // get n_rows and n_cols
1033 omb.n_rows = calc_map_dims(min_lat, max_lat);
1034 omb.n_cols = calc_map_dims(min_lon, max_lon);
1035 }
1036
1037 else {
1038 // Ensure odd dimensions
1039 omb.n_rows = (omb.wcs.naxis[1] % 2 == 0) ? omb.wcs.naxis[1] + 1 : omb.wcs.naxis[1];
1040 omb.n_cols = (omb.wcs.naxis[0] % 2 == 0) ? omb.wcs.naxis[0] + 1 : omb.wcs.naxis[0];
1041 }
1042
1043 Eigen::VectorXd rows_tan_vec = Eigen::VectorXd::LinSpaced(omb.n_rows, 0, omb.n_rows - 1).array() * omb.pixel_size_rad -
1044 (omb.n_rows / 2.0) * omb.pixel_size_rad;
1045 Eigen::VectorXd cols_tan_vec = Eigen::VectorXd::LinSpaced(omb.n_cols, 0, omb.n_cols - 1).array() * omb.pixel_size_rad -
1046 (omb.n_cols / 2.0) * omb.pixel_size_rad;
1047
1048
1049 // push back map sizes and coordinates
1050 map_extents.push_back({static_cast<int>(omb.n_rows), static_cast<int>(omb.n_cols)});
1051 map_coords.push_back({std::move(rows_tan_vec), std::move(cols_tan_vec)});
1052}
1053
1054// determine the map dimensions of the coadded map buffer
1055template <class EngineType>
1056void TimeOrderedDataProc<EngineType>::calc_cmb_size(std::vector<map_coord_t> &map_coords) {
1057 auto& cmb = engine().cmb;
1058
1059 // Initialize min/max values
1060 double min_row = std::numeric_limits<double>::max();
1061 double max_row = std::numeric_limits<double>::lowest();
1062 double min_col = min_row;
1063 double max_col = max_row;
1064
1065 // Find global min/max for rows and columns
1066 for (const auto& coord : map_coords) {
1067 min_row = std::min(min_row, coord.front().minCoeff());
1068 max_row = std::max(max_row, coord.front().maxCoeff());
1069 min_col = std::min(min_col, coord.back().minCoeff());
1070 max_col = std::max(max_col, coord.back().maxCoeff());
1071 }
1072
1073 // calculate dimensions
1074 auto calc_map_dims = [&](auto min_dim, auto max_dim) {
1075 int min_pix = static_cast<int>(ceil(abs(min_dim / engine().cmb.pixel_size_rad)));
1076 int max_pix = static_cast<int>(ceil(abs(max_dim / engine().cmb.pixel_size_rad)));
1077
1078 int n_dim = 2 * std::max(min_pix, max_pix) + 1;
1079 Eigen::VectorXd dim_vec = Eigen::VectorXd::LinSpaced(n_dim, 0, n_dim - 1)
1080 .array() * engine().cmb.pixel_size_rad - (n_dim / 2.0) * engine().cmb.pixel_size_rad;
1081
1082 return std::make_tuple(n_dim, std::move(dim_vec));
1083 };
1084
1085 // get dimensions and tangent coordinate vectorx
1086 auto [n_rows, rows_tan_vec] = calc_map_dims(min_row, max_row);
1087 auto [n_cols, cols_tan_vec] = calc_map_dims(min_col, max_col);
1088
1089 // Set dimensions and wcs parameters
1090 cmb.n_rows = n_rows;
1091 cmb.n_cols = n_cols;
1092 cmb.wcs.naxis[0] = n_cols;
1093 cmb.wcs.naxis[1] = n_rows;
1094 cmb.wcs.crpix[0] = (n_cols - 1) / 2.0;
1095 cmb.wcs.crpix[1] = (n_rows - 1) / 2.0;
1096 // set tangent plane coordinate vectors
1097 cmb.rows_tan_vec = std::move(rows_tan_vec);
1098 cmb.cols_tan_vec = std::move(cols_tan_vec);
1099}
1100
1101// allocate observation map buffer
1102template <class EngineType>
1104 auto& omb = engine().omb;
1105
1106 std::vector<Eigen::MatrixXd>().swap(omb.signal);
1107 std::vector<Eigen::MatrixXd>().swap(omb.weight);
1108 std::vector<Eigen::MatrixXd>().swap(omb.kernel);
1109 std::vector<Eigen::MatrixXd>().swap(omb.coverage);
1110 std::vector<Eigen::MatrixXd>().swap(omb.pointing);
1111
1112 // set omb dimensions and wcs parameters
1113 omb.n_rows = map_extent[0];
1114 omb.n_cols = map_extent[1];
1115 omb.wcs.naxis[0] = omb.n_cols;
1116 omb.wcs.naxis[1] = omb.n_rows;
1117 omb.wcs.crpix[0] = (omb.n_cols - 1) / 2.0;
1118 omb.wcs.crpix[1] = (omb.n_rows - 1) / 2.0;
1119 // set tangent plane coordinate vectors
1120 omb.rows_tan_vec = map_coord[0];
1121 omb.cols_tan_vec = map_coord[1];
1122
1123 Eigen::MatrixXd zero_matrix = Eigen::MatrixXd::Zero(omb.n_rows, omb.n_cols);
1124
1125 for (Eigen::Index i=0; i<engine().n_maps; ++i) {
1126 omb.signal.push_back(zero_matrix);
1127 omb.weight.push_back(zero_matrix);
1128
1129 if (engine().rtcproc.run_kernel) {
1130 omb.kernel.push_back(zero_matrix);
1131 }
1132
1133 if (engine().map_grouping != "detector") {
1134 omb.coverage.push_back(zero_matrix);
1135 }
1136 }
1137
1138 if (engine().rtcproc.run_polarization) {
1139 // allocate pointing matrix
1140 for (Eigen::Index i=0; i<engine().n_maps/engine().rtcproc.polarization.stokes_params.size(); ++i) {
1141 omb.pointing.emplace_back(omb.n_rows*omb.n_cols, 9);
1142 engine().omb.pointing.back().setZero();
1143 }
1144 }
1145}
1146
1147// allocate the coadded map buffer
1148template <class EngineType>
1150 auto& cmb = engine().cmb;
1151
1152 // clear map vectors
1153 std::vector<Eigen::MatrixXd>().swap(cmb.signal);
1154 std::vector<Eigen::MatrixXd>().swap(cmb.weight);
1155 std::vector<Eigen::MatrixXd>().swap(cmb.kernel);
1156 std::vector<Eigen::MatrixXd>().swap(cmb.coverage);
1157 std::vector<Eigen::MatrixXd>().swap(cmb.pointing);
1158
1159 Eigen::MatrixXd zero_matrix = Eigen::MatrixXd::Zero(cmb.n_rows, cmb.n_cols);
1160
1161 // loop through maps and allocate space
1162 for (Eigen::Index i=0; i<engine().n_maps; ++i) {
1163 cmb.signal.push_back(zero_matrix);
1164 cmb.weight.push_back(zero_matrix);
1165
1166 if (engine().rtcproc.run_kernel) {
1167 // allocate kernel
1168 cmb.kernel.push_back(zero_matrix);
1169 }
1170 if (engine().map_grouping!="detector") {
1171 // allocate coverage
1172 cmb.coverage.push_back(zero_matrix);
1173 }
1174 }
1175
1176 if (engine().rtcproc.run_polarization) {// && engine().run_noise) {
1177 // allocate pointing matrix
1178 for (Eigen::Index i=0; i<engine().n_maps/engine().rtcproc.polarization.stokes_params.size(); ++i) {
1179 cmb.pointing.emplace_back(cmb.n_rows*cmb.n_cols, 9);
1180 cmb.pointing.back().setZero();
1181 }
1182 }
1183}
1184
1185template <class EngineType>
1186template <class map_buffer_t>
1188 // clear noise map buffer
1189 std::vector<Eigen::Tensor<double,3>>().swap(nmb.noise);
1190
1191 // resize noise maps (n_maps, [n_rows, n_cols, n_noise])
1192 for (Eigen::Index i=0; i<engine().n_maps; ++i) {
1193 nmb.noise.emplace_back(nmb.n_rows, nmb.n_cols, nmb.n_noise);
1194 nmb.noise.at(i).setZero();
1195 }
1196}
1197
1198// coadd maps
1199template <class EngineType>
1201 // calculate the offset between cmb and omb
1202 int delta_row = 0.5*(engine().cmb.n_rows - engine().omb.n_rows);
1203 int delta_col= 0.5*(engine().cmb.n_cols - engine().omb.n_cols);
1204
1205 // loop through the maps
1206 for (Eigen::Index i=0; i<engine().n_maps; ++i) {
1207 // define common block references
1208 auto cmb_weight_block = engine().cmb.weight.at(i).block(delta_row, delta_col, engine().omb.n_rows, engine().omb.n_cols);
1209 auto cmb_signal_block = engine().cmb.signal.at(i).block(delta_row, delta_col, engine().omb.n_rows, engine().omb.n_cols);
1210
1211 // update cmb.weight with omb.weight
1212 cmb_weight_block += engine().omb.weight.at(i);
1213
1214 // update cmb.signal with omb.signal * omb.weight
1215 cmb_signal_block += (engine().omb.signal.at(i).array() * engine().omb.weight.at(i).array()).matrix();
1216
1217 // update cmb.kernel with omb.kernel * omb.weight
1218 if (engine().rtcproc.run_kernel) {
1219 auto cmb_kernel_block = engine().cmb.kernel.at(i).block(delta_row, delta_col, engine().omb.n_rows, engine().omb.n_cols);
1220 cmb_kernel_block += (engine().omb.kernel.at(i).array() * engine().omb.weight.at(i).array()).matrix();
1221 }
1222
1223 // update coverage
1224 if (!engine().cmb.coverage.empty()) {
1225 auto cmb_coverage_block = engine().cmb.coverage.at(i).block(delta_row, delta_col, engine().omb.n_rows, engine().omb.n_cols);
1226 cmb_coverage_block += engine().omb.coverage.at(i);
1227 }
1228 }
1229}
1230
1231template <class EngineType>
1233 // clear fits_io vectors
1234 std::vector<fitsIO<file_type_enum::write_fits, CCfits::ExtHDU*>>().swap(engine().coadd_fits_io_vec);
1235 std::vector<fitsIO<file_type_enum::write_fits, CCfits::ExtHDU*>>().swap(engine().coadd_noise_fits_io_vec);
1236 std::vector<fitsIO<file_type_enum::write_fits, CCfits::ExtHDU*>>().swap(engine().filtered_coadd_fits_io_vec);
1237 std::vector<fitsIO<file_type_enum::write_fits, CCfits::ExtHDU*>>().swap(engine().filtered_coadd_noise_fits_io_vec);
1238
1239 // loop through arrays
1240 for (Eigen::Index i=0; i<engine().calib.n_arrays; ++i) {
1241 // array index
1242 auto array = engine().calib.arrays[i];
1243 // array name
1244 std::string array_name = engine().toltec_io.array_name_map[array];
1245 // map filename
1246 auto filename = engine().toltec_io.template create_filename<engine_utils::toltecIO::toltec, engine_utils::toltecIO::map,
1247 engine_utils::toltecIO::raw>(engine().coadd_dir_name + "raw/",
1248 "", array_name, "",
1249 engine().telescope.sim_obs);
1250 // create fits_io class for current array file
1252 // append to fits_io vector
1253 engine().coadd_fits_io_vec.push_back(std::move(fits_io));
1254
1255 // if noise maps requested
1256 if (engine().run_noise) {
1257 // noise map filename
1258 auto filename = engine().toltec_io.template create_filename<engine_utils::toltecIO::toltec, engine_utils::toltecIO::noise,
1259 engine_utils::toltecIO::raw>(engine().coadd_dir_name + "raw/",
1260 "", array_name, "",
1261 engine().telescope.sim_obs);
1262 // create fits_io class for current array file
1264 // append to fits_io vector
1265 engine().coadd_noise_fits_io_vec.push_back(std::move(fits_io));
1266 }
1267 }
1268
1269 // if map filtering are requested
1270 if (engine().run_map_filter) {
1271 // loop through arrays
1272 for (Eigen::Index i=0; i<engine().calib.n_arrays; ++i) {
1273 // array index
1274 auto array = engine().calib.arrays[i];
1275 // array name
1276 std::string array_name = engine().toltec_io.array_name_map[array];
1277 // filtered map filename
1278 auto filename = engine().toltec_io.template create_filename<engine_utils::toltecIO::toltec, engine_utils::toltecIO::map,
1279 engine_utils::toltecIO::filtered>(engine().coadd_dir_name +
1280 "filtered/","", array_name,
1281 "", engine().telescope.sim_obs);
1282 // create fits_io class for current array file
1284 // append to fits_io vector
1285 engine().filtered_coadd_fits_io_vec.push_back(std::move(fits_io));
1286
1287 // if noise maps requested
1288 if (engine().run_noise) {
1289 // filtered noise map filename
1290 auto filename = engine().toltec_io.template create_filename<engine_utils::toltecIO::toltec, engine_utils::toltecIO::noise,
1291 engine_utils::toltecIO::filtered>(engine().coadd_dir_name +
1292 "filtered/","", array_name,
1293 "", engine().telescope.sim_obs);
1294 // create fits_io class for current array file
1296 // append to fits_io vector
1297 engine().filtered_coadd_noise_fits_io_vec.push_back(std::move(fits_io));
1298 }
1299 }
1300 }
1301}
1302
1303template <class EngineType>
1305 // get sortedfiles and directories in filepath
1306 std::set<fs::path> sorted_by_name;
1307 for (auto &entry : fs::directory_iterator(filepath))
1308 sorted_by_name.insert(entry);
1309
1310 // yaml node to store names
1311 YAML::Node node;
1312 // data products
1313 node["description"].push_back("citlali data products");
1314 // datetime when file is created
1315 node["date"].push_back(engine_utils::current_date_time());
1316 // citlali version
1317 node["citlali_version"].push_back(CITLALI_GIT_VERSION);
1318 // kids version
1319 node["kids_version"].push_back(KIDSCPP_GIT_VERSION);
1320 // tula version
1321 node["tula_version"].push_back(TULA_GIT_VERSION);
1322
1323 // call make_index_file recursively if current object is directory
1324 for (const auto & entry : sorted_by_name) {
1325 std::string path_string{entry.generic_string()};
1326 if (fs::is_directory(entry)) {
1327 make_index_file(path_string);
1328 }
1329 node["files/dirs"].push_back(path_string.substr(path_string.find_last_of("/") + 1));
1330 }
1331 // output yaml index file
1332 std::ofstream fout(filepath + "/index.yaml");
1333 fout << node;
1334}
Definition engine.h:153
@ map
Definition toltec_io.h:22
@ noise
Definition toltec_io.h:23
@ toltec
Definition toltec_io.h:14
@ filtered
Definition toltec_io.h:34
@ raw
Definition toltec_io.h:33
Definition fits_io.h:10
config::ConfigValidatorMixin< Derived, config::YamlConfig > ConfigMapper
Definition main_old.cpp:124
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
static long long modified_julian_date_to_unix(double jd)
Definition utils.h:157
static const std::string current_date_time()
Definition utils.h:91
Definition calib.h:12
Definition timestream.h:59
Definition netcdf_io.h:5
Definition todproc.h:16
friend OStream & operator<<(OStream &os, const DummyEngine &e)
Definition todproc.h:18
The DataItem struct This represent a single data item that belongs to a particular observation.
Definition io.h:58
The raw obs struct This represents a single observation that contains a set of data items and calibra...
Definition io.h:50
auto kidsdata() const -> decltype(auto)
Definition io.h:283
The time ordered data processing struct This wraps around the lali config.
Definition todproc.h:29
void calc_map_num()
Definition todproc.h:797
Eigen::MatrixXI scanindicies_t
Definition todproc.h:33
const Engine & engine() const
Definition todproc.h:126
void get_apt_from_files(const RawObs &rawobs)
Definition todproc.h:141
void make_index_file(std::string)
Definition todproc.h:1304
void get_tone_freqs_from_files(const RawObs &rawobs)
Definition todproc.h:217
void calc_omb_size(std::vector< map_extent_t > &, std::vector< map_coord_t > &)
Definition todproc.h:932
ConfigMapper< TimeOrderedDataProc< EngineType > > Base
Definition todproc.h:30
void allocate_omb(map_extent_t &, map_coord_t &)
Definition todproc.h:1103
void check_inputs(const RawObs &rawobs)
Definition todproc.h:406
void create_coadded_map_files()
Definition todproc.h:1232
void coadd()
Definition todproc.h:1200
void create_output_dir()
Definition todproc.h:347
std::size_t map_count_t
Definition todproc.h:36
void interp_pointing()
Definition todproc.h:737
void calc_cmb_size(std::vector< map_coord_t > &)
Definition todproc.h:1056
std::vector< std::tuple< Eigen::Index, Eigen::Index > > array_indices_t
Definition todproc.h:37
std::vector< Eigen::VectorXd > map_coord_t
Definition todproc.h:35
Engine m_engine
Definition todproc.h:136
TimeOrderedDataProc(config_t config)
Definition todproc.h:43
void align_timestreams_2(const RawObs &rawobs)
Definition todproc.h:734
std::shared_ptr< spdlog::logger > logger
Definition todproc.h:41
static auto check_config(const config_t &config) -> std::optional< std::string >
Definition todproc.h:46
Engine & engine()
Definition todproc.h:124
typename Base::config_t config_t
Definition todproc.h:31
std::vector< std::tuple< Eigen::Index, Eigen::Index > > det_indices_t
Definition todproc.h:38
void allocate_cmb()
Definition todproc.h:1149
void get_adc_snap_from_files(const RawObs &rawobs)
Definition todproc.h:311
std::vector< int > map_extent_t
Definition todproc.h:34
friend OStream & operator<<(OStream &os, const TimeOrderedDataProc &todproc)
Definition todproc.h:129
void allocate_nmb(map_buffer_t &)
Definition todproc.h:1187
void align_timestreams(const RawObs &rawobs)
Definition todproc.h:439