Citlali
Loading...
Searching...
No Matches
kernel.h
Go to the documentation of this file.
1#pragma once
2
3#include <string>
4
5#include <boost/math/special_functions/bessel.hpp>
6
10
11namespace timestream {
12
13class Kernel {
14public:
15 // get logger
16 std::shared_ptr<spdlog::logger> logger = spdlog::get("citlali_logger");
17
18 std::string filepath, type;
19 std::vector<std::string> img_ext_names;
20
21 // input kernel images
22 std::vector<Eigen::MatrixXd> images;
23
24 // sigma and fwhm from config
26
27 // limit on distance to calc sigma to
28 double sigma_limit = 3;
29
30 // map grouping
31 std::string map_grouping;
32
33 // initial setup
34 void setup(Eigen::Index, std::string);
35
36 // symmetric gaussian kernel
37 template<typename apt_t>
39 apt_t &);
40 // asymmetric elliptical gaussian kernel
41 template<typename apt_t>
43 // airy pattern kernel
44 template<typename apt_t>
46
47 // kernel from fits file
48 template<typename apt_t, typename Derived>
50 apt_t &, double, Eigen::DenseBase<Derived> &);
51};
52
53void Kernel::setup(Eigen::Index n_maps, std::string pixel_axes) {
54 if (type == "fits") {
55 if (img_ext_names.size()!=n_maps && img_ext_names.size()!=1) {
56 logger->error("mismatch for number of kernel images");
57 std::exit(EXIT_FAILURE);
58 }
59
61 for (auto & img_ext_name : img_ext_names) {
62 images.push_back(fits_io.get_hdu(img_ext_name));
63 }
64 }
65}
66
67template<typename apt_t>
69
70 // dimensions of scan
71 Eigen::Index n_dets = in.scans.data.cols();
72 Eigen::Index n_pts = in.scans.data.rows();
73
74 // resize kernel to match data size
75 in.kernel.data.resize(n_pts,n_dets);
76
77 double sigma = sigma_rad;
78
79 for (Eigen::Index i=0; i<n_dets; ++i) {
80 // calc tangent plane pointing
81 auto [lat, lon] = engine_utils::calc_det_pointing(in.tel_data.data, apt["x_t"](i), apt["y_t"](i),
82 pixel_axes, in.pointing_offsets_arcsec.data, map_grouping);
83
84 // distance to source to truncate it
85 auto dist = ((lat.array()).pow(2) + (lon.array()).pow(2)).sqrt();
86
87 // calculate stddev from apt table if config stddev <=0
88 if (sigma_rad <= 0) {
89 sigma = FWHM_TO_STD * ASEC_TO_RAD*(apt["a_fwhm"](i) + apt["b_fwhm"](i))/2;
90 }
91
92 // loop through samples and calculate
93 for (Eigen::Index j=0; j<n_pts; ++j) {
94 // truncate within radius
95 if (dist(j) <= sigma_limit*sigma) {
96 in.kernel.data(j,i) = exp(-0.5*pow(dist(j)/sigma,2));
97 }
98 else {
99 in.kernel.data(j,i) = 0;
100 }
101 }
102 }
103}
104
105template<typename apt_t>
106void Kernel::create_gaussian_kernel(TCData<TCDataKind::RTC, Eigen::MatrixXd> &in, std::string &pixel_axes, apt_t &apt) {
107
108 // dimensions of scan
109 Eigen::Index n_dets = in.scans.data.cols();
110 Eigen::Index n_pts = in.scans.data.rows();
111
112 // resize kernel to match data size
113 in.kernel.data.resize(n_pts,n_dets);
114
115 // get parameters for current detector
116 double amp = 1.0;
117 double off_lat = 0.0;
118 double off_lon = 0.0;
119
120 // beam standard deviations
121 double sigma_lat = sigma_rad;
122 double sigma_lon = sigma_rad;
123
124 for (Eigen::Index i=0; i<n_dets; ++i) {
125 // calc tangent plane pointing
126 auto [lat, lon] = engine_utils::calc_det_pointing(in.tel_data.data, apt["x_t"](i), apt["y_t"](i),
127 pixel_axes, in.pointing_offsets_arcsec.data, map_grouping);
128
129 // distance to source to truncate it
130 auto dist = ((lat.array()).pow(2) + (lon.array()).pow(2)).sqrt();
131
132 // calculate stddev from apt table if config stddev <=0
133 if (sigma_rad <= 0) {
134 sigma_lat = FWHM_TO_STD * ASEC_TO_RAD * apt["b_fwhm"](i);
135 sigma_lon = FWHM_TO_STD * ASEC_TO_RAD * apt["a_fwhm"](i);
136 }
137
138 // rotation angle
139 double rot_ang = apt["angle"](i);
140
141 auto cost2 = cos(rot_ang) * cos(rot_ang);
142 auto sint2 = sin(rot_ang) * sin(rot_ang);
143 auto sin2t = sin(2. * rot_ang);
144 auto xstd2 = sigma_lon * sigma_lon;
145 auto ystd2 = sigma_lat * sigma_lat;
146 auto a = - 0.5 * ((cost2 / xstd2) + (sint2 / ystd2));
147 auto b = - 0.5 * ((sin2t / xstd2) - (sin2t / ystd2));
148 auto c = - 0.5 * ((sint2 / xstd2) + (cost2 / ystd2));
149
150 double sigma_limit_det = sigma_limit * (sigma_lat + sigma_lon)/2;
151
152 // make elliptical gaussian
153 for (Eigen::Index j=0; j<n_pts; ++j) {
154 // truncate within radius
155 if (dist(j) <= sigma_limit_det) {
156 in.kernel.data(j,i) = amp*exp(pow(lon(j) - off_lon, 2) * a +
157 (lon(j) - off_lon) * (lat(j) - off_lat) * b +
158 pow(lat(j) - off_lat, 2) * c);
159 }
160 else {
161 in.kernel.data(j,i) = 0;
162 }
163 }
164 }
165}
166
167template<typename apt_t>
168void Kernel::create_airy_kernel(TCData<TCDataKind::RTC, Eigen::MatrixXd> &in, std::string &pixel_axes, apt_t &apt) {
169
170 Eigen::Index n_dets = in.scans.data.cols();
171 Eigen::Index n_pts = in.scans.data.rows();
172
173 in.kernel.data.resize(n_pts,n_dets);
174
175 double fwhm = fwhm_rad;
176
177 // loop through detectors
178 for (Eigen::Index i=0; i<n_dets; ++i) {
179 // calc tangent plane pointing
180 auto [lat, lon] = engine_utils::calc_det_pointing(in.tel_data.data, apt["x_t"](i), apt["y_t"](i),
181 pixel_axes, in.pointing_offsets_arcsec.data, map_grouping);
182
183 // distance to source to truncate it
184 auto dist = ((lat.array()).pow(2) + (lon.array()).pow(2)).sqrt();
185
186 // get fwhm from apt if config file fwhm is <= 0
187 if (fwhm_rad <= 0) {
188 fwhm = ASEC_TO_RAD*(apt["a_fwhm"](i) + apt["b_fwhm"](i))/2;
189 }
190
191 // airy pattern factor
192 double factor = pi*(1.028/fwhm);
193
194 for (Eigen::Index j=0; j<n_pts; ++j) {
195 if (dist(j) <= sigma_limit*fwhm) {
196 in.kernel.data(j,i) = pow(2*boost::math::cyl_bessel_j(1,factor*dist(j))/(factor*dist(j)),2);
197 }
198 else {
199 in.kernel.data(j,i) = 0;
200 }
201 }
202 }
203}
204
205template<typename apt_t, typename Derived>
207 double pixel_size_rad, Eigen::DenseBase<Derived> &map_indices) {
208
209 Eigen::Index n_dets = in.scans.data.cols();
210 Eigen::Index n_pts = in.scans.data.rows();
211
212 in.kernel.data.resize(n_pts,n_dets);
213
214 Eigen::Index map_index = 0;
215
216 // loop through detectors
217 for (Eigen::Index i=0; i<n_dets; ++i) {
218 // calc tangent plane pointing
219 auto [lat, lon] = engine_utils::calc_det_pointing(in.tel_data.data, apt["x_t"](i), apt["y_t"](i),
220 pixel_axes, in.pointing_offsets_arcsec.data, map_grouping);
221
222 if (images.size() > 1) {
223 map_index = map_indices(i);
224 }
225
226 // get map buffer row and col indices for lat and lon vectors
227 Eigen::VectorXd irows = lat.array()/pixel_size_rad + (images[map_index].rows())/2.;
228 Eigen::VectorXd icols = lon.array()/pixel_size_rad + (images[map_index].cols())/2.;
229
230 for (Eigen::Index j = 0; j<n_pts; ++j) {
231 // row and col pixel for kernel image
232 Eigen::Index ir = irows(j);
233 Eigen::Index ic = icols(j);
234
235 // check if current sample is on the image and add to the timestream
236 if ((ir >= 0) && (ir < images[map_index].rows()) && (ic >= 0) && (ic < images[map_index].cols())) {
237 in.kernel.data(j,i) = images[map_index](ir,ic);
238 }
239 }
240 }
241}
242} // namespace timestream
Definition fits_io.h:10
auto get_hdu(std::string hdu_name)
Definition fits_io.h:90
Definition kernel.h:13
void setup(Eigen::Index, std::string)
Definition kernel.h:53
std::string type
Definition kernel.h:18
double sigma_rad
Definition kernel.h:25
double sigma_limit
Definition kernel.h:28
std::string filepath
Definition kernel.h:18
std::vector< std::string > img_ext_names
Definition kernel.h:19
void create_kernel_from_fits(TCData< TCDataKind::RTC, Eigen::MatrixXd > &, std::string &, apt_t &, double, Eigen::DenseBase< Derived > &)
Definition kernel.h:206
void create_symmetric_gaussian_kernel(TCData< TCDataKind::RTC, Eigen::MatrixXd > &, std::string &, apt_t &)
Definition kernel.h:68
std::vector< Eigen::MatrixXd > images
Definition kernel.h:22
std::string map_grouping
Definition kernel.h:31
void create_gaussian_kernel(TCData< TCDataKind::RTC, Eigen::MatrixXd > &, std::string &, apt_t &)
Definition kernel.h:106
double fwhm_rad
Definition kernel.h:25
std::shared_ptr< spdlog::logger > logger
Definition kernel.h:16
void create_airy_kernel(TCData< TCDataKind::RTC, Eigen::MatrixXd > &, std::string &, apt_t &)
Definition kernel.h:168
#define FWHM_TO_STD
Definition constants.h:48
#define ASEC_TO_RAD
Definition constants.h:33
constexpr auto pi
Definition constants.h:6
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
Definition clean.h:10
TC data class.
Definition timestream.h:55