Citlali
Loading...
Searching...
No Matches
calibrate.h
Go to the documentation of this file.
1#pragma once
2
3#include <map>
4#include <string>
5#include <Eigen/Core>
6
9
10namespace timestream {
11
13public:
14 // get logger
15 std::shared_ptr<spdlog::logger> logger = spdlog::get("citlali_logger");
16
17 // extinction model
18 std::string extinction_model;
19
20 // coefficients for transmission ratios fit to order 6 polynomial
21 std::map<std::string,Eigen::VectorXd> tx_ratio_coeff;
22
23 std::map<std::string, double> tx_225_zenith = {
24 {"am_q25",0.9500275},
25 {"am_q50",0.9142065},
26 {"am_q75",0.8515054},
27 {"am_q95",0.7337698},
28 };
29
30 void setup(double tau_225_zenith) {
31 // cos of zenith angle
32 auto cz = cos(pi/2 - 80.0*DEG_TO_RAD);
33 // 1/cos(zenith angle)
34 auto secz = 1. / cz;
35 // airmass
36 auto A = secz * (1. - 0.0012 * (pow(secz, 2) - 1.));
37
38 // tau at 225 GHz at 80deg from atm model and airmass
39 Eigen::VectorXd tau_225_calc(tx_225_zenith.size());
40
41 // calc tau at 225 GHz at 80deg from each model
42 int i = 0;
43 for (const auto &[key,val]: tx_225_zenith) {
44 tau_225_calc(i) = -log(val)/A;
45 i++;
46 }
47
48 // set initial model to am_q25
49 extinction_model = "am_q25";
50
51 // find model with closest tau to telescope tau and use that model
52 // for extinction correction
53 i = 0;
54 for (const auto &[key,val]: tx_225_zenith) {
55 if (tau_225_calc(i) <= tau_225_zenith) {
56 extinction_model = key;
57 }
58 i++;
59 }
60
61 // allocate transmission coefficients (order 6 polynomial)
62 tx_ratio_coeff["a1100"].resize(7);
63 tx_ratio_coeff["a1400"].resize(7);
64 tx_ratio_coeff["a2000"].resize(7);
65
66 // am_q25
67 if (extinction_model=="am_q25") {
68 tx_ratio_coeff["a1100"] << -0.12008024, 0.72422015, -1.81734478, 2.45313012, -1.92159695, 0.86918801, 0.78604295;
69 tx_ratio_coeff["a1400"] << 0.02619509, -0.15757661, 0.39400473, -0.52912696, 0.411213, -0.18360141, 1.04398466;
70 tx_ratio_coeff["a2000"] << 0.16726241, -1.00436302, 2.50507317, -3.35219659, 2.59080373, -1.14622096, 1.26931683;
71 }
72
73 // am_q50
74 else if (extinction_model=="am_q50") {
75 tx_ratio_coeff["a1100"] << -0.18770822, 1.13390437, -2.85173457, 3.8617083, -3.03996805, 1.38624303,
76 0.65300169;
77 tx_ratio_coeff["a1400"] << 0.04292884, -0.25817762, 0.64533115, -0.86622214, 0.67267823, -0.29996916,
78 1.07167603;
79 tx_ratio_coeff["a2000"] << 0.35178447, -2.10859714, 5.24620825, -6.9952531, 5.37645792, -2.35675076,
80 1.54286813;
81 }
82
83 // am_q75
84 else if (extinction_model=="am_q75") {
85 tx_ratio_coeff["a1100"] << -0.28189529, 1.70842347, -4.31606883, 5.88248549, -4.67702093, 2.16747228,
86 0.4393435;
87 tx_ratio_coeff["a1400"] << 0.07581154, -0.45574885, 1.13852458, -1.52697451, 1.18425865, -0.5269455,
88 1.12531753;
89 tx_ratio_coeff["a2000"] << 0.79869908, -4.77265095, 11.82393401, -15.67007557, 11.93052031,
90 -5.14788907, 2.14595898;
91 }
92
93 // am_q95
94 else if (extinction_model=="am_q95") {
95 tx_ratio_coeff["a1100"] << -1.21882233, 6.67068453, -14.96466875, 17.78045563, -12.10288687,
96 4.76050807, -0.06765066;
97 tx_ratio_coeff["a1400"] << 0.76090502, -4.05867663, 8.78487281, -9.90872343, 6.2198602, -2.13790165,
98 1.3668983;
99 tx_ratio_coeff["a2000"] << 16.0063036, -84.30325144, 179.28096414, -197.05751682, 118.73627425,
100 -37.99279818, 6.55457576;
101 }
102 }
103
104 // polynomial fit to transmission ratio as a function of elevation (radians)
105 template <typename DerivedA, typename DerivedB>
106 auto tau_polynomial(Eigen::DenseBase<DerivedA> &coeff, Eigen::DenseBase<DerivedB> &elev) {
107
108 return coeff(0)*pow(elev.derived().array(),6) + coeff(1)*pow(elev.derived().array(),5) +
109 coeff(2)*pow(elev.derived().array(),4) + coeff(3)*pow(elev.derived().array(),3) +
110 coeff(4)*pow(elev.derived().array(),2) + coeff(5)*pow(elev.derived().array(),1) +
111 coeff(6);
112 }
113
114 // calculate tau for each toltec band for given elevations
115 template <typename Derived>
116 auto calc_tau(Eigen::DenseBase<Derived> &, double);
117
118 // run flux calibration on the timestreams
119 template <TCDataKind tcdata_kind, class calib_t>
121
122 // run extinction correction on the timestreams
123 template <TCDataKind tcdata_kind, class calib_t, typename tau_t>
125};
126
127template <typename Derived>
128auto Calibration::calc_tau(Eigen::DenseBase<Derived> &elev, double tau_225_GHz) {
129
130 // tau at toltec freqs
131 std::map<int,Eigen::VectorXd> tau_freq;
132
133 // zenith angle
134 auto cz = cos(pi/2 - elev.derived().array());
135 auto secz = 1. / cz;
136 auto A = secz * (1. - 0.0012 * (pow(secz.array(), 2) - 1.));
137
138 // observed 225 GHz tau
139 auto obs_tau_i = A*tau_225_GHz;
140
141 // 225 GHz transmission
142 auto tx = (-obs_tau_i).exp();
143
144 // a1100
145 auto tx_a1100 = tau_polynomial(tx_ratio_coeff["a1100"],elev).array()*tx.array();
146 tau_freq[0] = -(tx_a1100.array().log())/A.array();
147
148 // a1400
149 auto tx_a1400 = tau_polynomial(tx_ratio_coeff["a1400"],elev).array()*tx.array();
150 tau_freq[1] = -(tx_a1400.array().log())/A.array();
151
152 // a2000
153 auto tx_a2000 = tau_polynomial(tx_ratio_coeff["a2000"],elev).array()*tx.array();
154 tau_freq[2] = -(tx_a2000.array().log())/A.array();
155
156 return tau_freq;
157}
158
159template <TCDataKind tcdata_kind, class calib_t>
161
162 // loop through detectors
163 for (Eigen::Index i=0; i<in.scans.data.cols(); ++i) {
164 // current array index in apt table
165 Eigen::Index array_index = calib.apt["array"](i);
166
167 // flux conversion factor for non-mJy/beam units
168 in.fcf.data(i) = calib.flux_conversion_factor(array_index);
169
170 // data x flxscale x factor
171 in.scans.data.col(i) = in.scans.data.col(i).array()*in.fcf.data(i)*calib.apt["flxscale"](i);
172 }
173}
174
175template <TCDataKind tcdata_kind, class calib_t, typename tau_t>
177
178 // loop through detectors
179 for (Eigen::Index i=0; i<in.scans.data.cols(); ++i) {
180 // current array index in apt table
181 Eigen::Index array_index = calib.apt["array"](i);
182
183 // factor = 1 / exp(-tau_freq)
184 auto factor = 1./(-tau_freq[array_index]).array().exp();
185
186 // apply extinction correction to fcf
187 in.fcf.data(i) = (in.fcf.data(i)*factor).mean();
188
189 // data x factor
190 in.scans.data.col(i) = in.scans.data.col(i).array()*factor.array();
191 }
192}
193} // namespace timestream
Definition calibrate.h:12
std::map< std::string, double > tx_225_zenith
Definition calibrate.h:23
void extinction_correction(TCData< tcdata_kind, Eigen::MatrixXd > &, calib_t &, tau_t &)
Definition calibrate.h:176
auto calc_tau(Eigen::DenseBase< Derived > &, double)
Definition calibrate.h:128
void setup(double tau_225_zenith)
Definition calibrate.h:30
std::shared_ptr< spdlog::logger > logger
Definition calibrate.h:15
std::map< std::string, Eigen::VectorXd > tx_ratio_coeff
Definition calibrate.h:21
auto tau_polynomial(Eigen::DenseBase< DerivedA > &coeff, Eigen::DenseBase< DerivedB > &elev)
Definition calibrate.h:106
void calibrate_tod(TCData< tcdata_kind, Eigen::MatrixXd > &, calib_t &)
Definition calibrate.h:160
std::string extinction_model
Definition calibrate.h:18
#define DEG_TO_RAD
Definition constants.h:27
constexpr auto pi
Definition constants.h:6
Definition clean.h:10
TC data class.
Definition timestream.h:55