tracking_validation.cpp 34 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589
  1. #include <iostream>
  2. #include <vector>
  3. #include <map>
  4. #include <utility>
  5. #include <numeric>
  6. #include <limits>
  7. #include <cmath>
  8. #include <TSystem.h>
  9. #include "filval/filval.hpp"
  10. #include "filval/root/filval.hpp"
  11. #include <boost/range/combine.hpp>
  12. #include <boost/format.hpp>
  13. #include "TrackingNtuple.h"
  14. #include "analysis/obj_types.cpp"
  15. using namespace std;
  16. using namespace fv;
  17. using namespace fv::root;
  18. typedef std::tuple<Seed, PixRecHit, Track, SimHit, SimTrack> MatchedTrack;
  19. typedef std::pair<std::vector<float>,std::vector<float>> pair_vec;
  20. #define HIT_TYPE_PIXEL 0
  21. #define HIT_TYPE_GLUED 1
  22. #define HIT_TYPE_STRIP 2
  23. #define PIXEL_BARREL 1
  24. #define PIXEL_ENDCAP 2
  25. std::map<int,string> subdet_names = {{1, "BPIX"}, {2, "FPIX"}};
  26. template<typename A, typename B>
  27. float displacement(const A& a, const B& b){
  28. return std::sqrt(pow(a.x-b.x, 2) +
  29. pow(a.x-b.x, 2) +
  30. pow(a.x-b.x, 2));
  31. }
  32. bool in_det(const MatchedTrack &mt, int &&det, int &&layer){
  33. auto& hit = std::get<PixRecHit>(mt);
  34. return hit.det == det && hit.lay == layer;
  35. };
  36. float pseudorapidity(float z, float r){
  37. float theta = atan2(r, z);
  38. return -log(tan(theta/2.0));
  39. }
  40. float pseudorapidity(float x, float y, float z){
  41. float r = sqrt(x*x + y*y);
  42. return pseudorapidity(z, r);
  43. }
  44. template<typename T>
  45. float pseudorapidity(const T &t){
  46. return pseudorapidity(t.x, t.y, t.z);
  47. }
  48. vector<string> seedTypes =
  49. {"initialStepSeeds",
  50. "highPtTripletStepSeeds",
  51. "mixedTripletStepSeeds",
  52. "pixelLessStepSeeds",
  53. "tripletElectronSeeds",
  54. "pixelPairElectronSeeds",
  55. "stripPairElectronSeeds"};
  56. Value<vector<Seed>>* seeds;
  57. Value<vector<PixRecHit>>* pixrec_hits;
  58. Value<vector<Track>>* tracks;
  59. Value<vector<SimHit>>* sim_hits;
  60. Value<vector<SimTrack>>* sim_tracks;
  61. void register_objects(TrackingDataSet& tds){
  62. seeds = register_seeds(tds);
  63. pixrec_hits = register_pixrec_hits(tds);
  64. tracks = register_tracks(tds);
  65. sim_hits = register_sim_hits(tds);
  66. sim_tracks = register_sim_tracks(tds);
  67. }
  68. void setup_first_hit_pairs(TrackingDataSet& tds){
  69. typedef std::tuple<PixRecHit, SimHit> HitPair;
  70. auto& find_matched_nth_hit_in_layer =
  71. func<vector<HitPair>(vector<Seed>,
  72. vector<PixRecHit>,
  73. vector<SimHit>,
  74. int, int, int)>("find_matched_nth_hit_in_layer",
  75. FUNC(([](const vector<Seed>& seeds,
  76. const vector<PixRecHit>& pixrec_hits,
  77. const vector<SimHit>& sim_hits,
  78. const int&& det,
  79. const int&& bpix_layer,
  80. const int&& hit_number){
  81. vector<HitPair> matched_hits;
  82. for(const Seed &seed : seeds){ // loop over all seeds
  83. if(seed.hitIdx.size() <= hit_number) continue;
  84. if(seed.algoOriginal < 0 || seed.algoOriginal >= seedTypes.size()) continue;
  85. if(seed.hitType[hit_number] != HIT_TYPE_PIXEL) continue; // take only pixel hits for now
  86. const PixRecHit &rec_hit = pixrec_hits[seed.hitIdx[hit_number]];
  87. if(rec_hit.det == det && rec_hit.lay == bpix_layer){
  88. if(rec_hit.simHitIdx.size() > 0){
  89. // take first matched simhit (should be the closest one)
  90. matched_hits.push_back({rec_hit, sim_hits[rec_hit.simHitIdx[0]]});
  91. }
  92. }
  93. }
  94. return matched_hits;
  95. })));
  96. auto& select_even_odd_ladder_blade_hit_pairs =
  97. func<vector<HitPair>(vector<HitPair>, bool)>("select_even_odd_ladder_blade_hit_pairs",
  98. FUNC(([](const vector<HitPair>& hit_pairs,
  99. const bool &&odd){
  100. vector<HitPair> even_pairs;
  101. for(const HitPair &hit_pair : hit_pairs){ // loop over all seeds
  102. if(std::get<PixRecHit>(hit_pair).ladder_blade % 2 == odd){
  103. even_pairs.push_back(hit_pair);
  104. }
  105. }
  106. return even_pairs;
  107. })));
  108. auto barrel_val = constant("PIXEL_BARREL", PIXEL_BARREL);
  109. auto endcap_val = constant("PIXEL_ENDCAP", PIXEL_ENDCAP);
  110. auto first_hit = constant("1st", 0);
  111. auto second_hit = constant("2nd", 1);
  112. // First hits on inner three bpix layers
  113. auto first_hits_in_B1 = fv::apply(find_matched_nth_hit_in_layer,
  114. fv::tuple(seeds, pixrec_hits, sim_hits, barrel_val, constant("L1", 1), first_hit), "first_hits_in_B1");
  115. auto first_hits_in_B2 = fv::apply(find_matched_nth_hit_in_layer,
  116. fv::tuple(seeds, pixrec_hits, sim_hits, barrel_val, constant("L2", 2), first_hit), "first_hits_in_B2");
  117. auto first_hits_in_B3 = fv::apply(find_matched_nth_hit_in_layer,
  118. fv::tuple(seeds, pixrec_hits, sim_hits, barrel_val, constant("L3", 3), first_hit), "first_hits_in_B3");
  119. // Second hits on outer three bpix layers
  120. auto second_hits_in_B2 = fv::apply(find_matched_nth_hit_in_layer,
  121. fv::tuple(seeds, pixrec_hits, sim_hits, barrel_val, constant("L2", 2), second_hit), "second_hits_in_B2");
  122. auto second_hits_in_B3 = fv::apply(find_matched_nth_hit_in_layer,
  123. fv::tuple(seeds, pixrec_hits, sim_hits, barrel_val, constant("L3", 3), second_hit), "second_hits_in_B3");
  124. auto second_hits_in_B4 = fv::apply(find_matched_nth_hit_in_layer,
  125. fv::tuple(seeds, pixrec_hits, sim_hits, barrel_val, constant("L4", 4), second_hit), "second_hits_in_B4");
  126. /* auto first_hits_in_F1 = fv::apply(find_matched_nth_hit_in_layer, */
  127. /* fv::tuple(seeds, pixrec_hits, sim_hits, endcap_val, constant("1", 1)), "first_hits_in_F1"); */
  128. /* auto first_hits_in_F2 = fv::apply(find_matched_nth_hit_in_layer, */
  129. /* fv::tuple(seeds, pixrec_hits, sim_hits, endcap_val, constant("2", 2)), "first_hits_in_F2"); */
  130. // Even vs Odd Ladders
  131. auto even = constant("even", false);
  132. auto odd = constant("odd", true);
  133. auto first_hits_in_B1_even_ladder = fv::apply(select_even_odd_ladder_blade_hit_pairs,
  134. fv::tuple(first_hits_in_B1, even), "first_hits_in_B1_even_ladder");
  135. auto first_hits_in_B1_odd_ladder = fv::apply(select_even_odd_ladder_blade_hit_pairs,
  136. fv::tuple(first_hits_in_B1, odd), "first_hits_in_B1_odd_ladder");
  137. auto first_hits_in_B2_even_ladder = fv::apply(select_even_odd_ladder_blade_hit_pairs,
  138. fv::tuple(first_hits_in_B2, even), "first_hits_in_B2_even_ladder");
  139. auto first_hits_in_B2_odd_ladder = fv::apply(select_even_odd_ladder_blade_hit_pairs,
  140. fv::tuple(first_hits_in_B2, odd), "first_hits_in_B2_odd_ladder");
  141. auto first_hits_in_B3_even_ladder = fv::apply(select_even_odd_ladder_blade_hit_pairs,
  142. fv::tuple(first_hits_in_B3, even), "first_hits_in_B3_even_ladder");
  143. auto first_hits_in_B3_odd_ladder = fv::apply(select_even_odd_ladder_blade_hit_pairs,
  144. fv::tuple(first_hits_in_B3, odd), "first_hits_in_B3_odd_ladder");
  145. //Plots for dPhi of collections defined above
  146. auto& calc_dphi_v_eta = func<pair_vec(vector<HitPair>)>("calc_dphi_v_eta",
  147. FUNC(([](const vector<HitPair>& hit_pairs){
  148. vector<float> dphis;
  149. vector<float> etas;
  150. for(auto hit_pair : hit_pairs){
  151. auto& pixrec_hit = std::get<PixRecHit>(hit_pair);
  152. auto& sim_hit = std::get<SimHit>(hit_pair);
  153. dphis.push_back(atan2(sim_hit.x, sim_hit.y) - atan2(pixrec_hit.x, pixrec_hit.y));
  154. etas.push_back(pseudorapidity(pixrec_hit));
  155. }
  156. return std::make_pair(etas, dphis);
  157. })));
  158. TH2Params params_dphi = {"$\\eta$", 100, -4, 4,
  159. "$\\Delta \\phi$(rad)", 75, -.0015, .0015};
  160. tds.register_container<ContainerTH2Many<float>>("dphi_v_eta_first_hits_in_B1",
  161. fv::apply(calc_dphi_v_eta, fv::tuple(first_hits_in_B1)),
  162. "First Hit in BPIX-L1", params_dphi);
  163. tds.register_container<ContainerTH2Many<float>>("dphi_v_eta_first_hits_in_B2",
  164. fv::apply(calc_dphi_v_eta, fv::tuple(first_hits_in_B2)),
  165. "First Hit in BPIX-L2", params_dphi);
  166. tds.register_container<ContainerTH2Many<float>>("dphi_v_eta_first_hits_in_B3",
  167. fv::apply(calc_dphi_v_eta, fv::tuple(first_hits_in_B3)),
  168. "First Hit in BPIX-L3", params_dphi);
  169. tds.register_container<ContainerTH2Many<float>>("dphi_v_eta_first_hits_in_B1_even_ladder",
  170. fv::apply(calc_dphi_v_eta, fv::tuple(first_hits_in_B1_even_ladder)),
  171. "First Hit in BPIX-L1 - Even Ladders", params_dphi);
  172. tds.register_container<ContainerTH2Many<float>>("dphi_v_eta_first_hits_in_B1_odd_ladder",
  173. fv::apply(calc_dphi_v_eta, fv::tuple(first_hits_in_B1_odd_ladder)),
  174. "First Hit in BPIX-L1 - Odd Ladders", params_dphi);
  175. tds.register_container<ContainerTH2Many<float>>("dphi_v_eta_first_hits_in_B2_even_ladder",
  176. fv::apply(calc_dphi_v_eta, fv::tuple(first_hits_in_B2_even_ladder)),
  177. "First Hit in BPIX-L2 - Even Ladders", params_dphi);
  178. tds.register_container<ContainerTH2Many<float>>("dphi_v_eta_first_hits_in_B2_odd_ladder",
  179. fv::apply(calc_dphi_v_eta, fv::tuple(first_hits_in_B2_odd_ladder)),
  180. "First Hit in BPIX-L2 - Odd Ladders", params_dphi);
  181. tds.register_container<ContainerTH2Many<float>>("dphi_v_eta_first_hits_in_B3_even_ladder",
  182. fv::apply(calc_dphi_v_eta, fv::tuple(first_hits_in_B3_even_ladder)),
  183. "First Hit in BPIX-L3 - Even Ladders", params_dphi);
  184. tds.register_container<ContainerTH2Many<float>>("dphi_v_eta_first_hits_in_B3_odd_ladder",
  185. fv::apply(calc_dphi_v_eta, fv::tuple(first_hits_in_B3_odd_ladder)),
  186. "First Hit in BPIX-L3 - Odd Ladders", params_dphi);
  187. tds.register_container<ContainerTH2Many<float>>("dphi_v_eta_second_hits_in_B2",
  188. fv::apply(calc_dphi_v_eta, fv::tuple(second_hits_in_B2)),
  189. "Second Hit in BPIX-L2", params_dphi);
  190. tds.register_container<ContainerTH2Many<float>>("dphi_v_eta_second_hits_in_B3",
  191. fv::apply(calc_dphi_v_eta, fv::tuple(second_hits_in_B3)),
  192. "Second Hit in BPIX-L3", params_dphi);
  193. tds.register_container<ContainerTH2Many<float>>("dphi_v_eta_second_hits_in_B4",
  194. fv::apply(calc_dphi_v_eta, fv::tuple(second_hits_in_B4)),
  195. "Second Hit in BPIX-L4", params_dphi);
  196. //Plots for dz of collections defined above
  197. auto& calc_dz_v_eta = func<pair_vec(vector<HitPair>)>("calc_dz_v_eta",
  198. FUNC(([](const vector<HitPair>& hit_pairs){
  199. vector<float> dzs;
  200. vector<float> etas;
  201. for(auto hit_pair : hit_pairs){
  202. auto& pixrec_hit = std::get<PixRecHit>(hit_pair);
  203. auto& sim_hit = std::get<SimHit>(hit_pair);
  204. dzs.push_back(sim_hit.z - pixrec_hit.z);
  205. etas.push_back(pseudorapidity(pixrec_hit));
  206. }
  207. return std::make_pair(etas, dzs);
  208. })));
  209. TH2Params params_dz = {"$\\eta$", 100, -4, 4,
  210. "$\\Delta z$(rad)", 100, -.01, .01};
  211. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_first_hits_in_B1",
  212. fv::apply(calc_dz_v_eta, fv::tuple(first_hits_in_B1)),
  213. "First Hit in BPIX-L1", params_dz);
  214. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_first_hits_in_B2",
  215. fv::apply(calc_dz_v_eta, fv::tuple(first_hits_in_B2)),
  216. "First Hit in BPIX-L2", params_dz);
  217. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_first_hits_in_B3",
  218. fv::apply(calc_dz_v_eta, fv::tuple(first_hits_in_B3)),
  219. "First Hit in BPIX-L3", params_dz);
  220. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_first_hits_in_B1_even_ladder",
  221. fv::apply(calc_dz_v_eta, fv::tuple(first_hits_in_B1_even_ladder)),
  222. "First Hit in BPIX-L1 - Even Ladders", params_dz);
  223. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_first_hits_in_B1_odd_ladder",
  224. fv::apply(calc_dz_v_eta, fv::tuple(first_hits_in_B1_odd_ladder)),
  225. "First Hit in BPIX-L1 - Odd Ladders", params_dz);
  226. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_first_hits_in_B2_even_ladder",
  227. fv::apply(calc_dz_v_eta, fv::tuple(first_hits_in_B2_even_ladder)),
  228. "First Hit in BPIX-L2 - Even Ladders", params_dz);
  229. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_first_hits_in_B2_odd_ladder",
  230. fv::apply(calc_dz_v_eta, fv::tuple(first_hits_in_B2_odd_ladder)),
  231. "First Hit in BPIX-L2 - Odd Ladders", params_dz);
  232. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_first_hits_in_B3_even_ladder",
  233. fv::apply(calc_dz_v_eta, fv::tuple(first_hits_in_B3_even_ladder)),
  234. "First Hit in BPIX-L3 - Even Ladders", params_dz);
  235. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_first_hits_in_B3_odd_ladder",
  236. fv::apply(calc_dz_v_eta, fv::tuple(first_hits_in_B3_odd_ladder)),
  237. "First Hit in BPIX-L3 - Odd Ladders", params_dz);
  238. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_second_hits_in_B2",
  239. fv::apply(calc_dz_v_eta, fv::tuple(second_hits_in_B2)),
  240. "Second Hit in BPIX-L2", params_dz);
  241. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_second_hits_in_B3",
  242. fv::apply(calc_dz_v_eta, fv::tuple(second_hits_in_B3)),
  243. "Second Hit in BPIX-L3", params_dz);
  244. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_second_hits_in_B4",
  245. fv::apply(calc_dz_v_eta, fv::tuple(second_hits_in_B4)),
  246. "Second Hit in BPIX-L4", params_dz);
  247. }
  248. void setup_matched_tracks(TrackingDataSet& tds){
  249. // Finds sets of (seed, rechit, track, simhit sim_track) where the cycle
  250. // track->seed->rec_hit->sim_hit->sim_track->track
  251. // is closed
  252. auto& find_matched_tracks =
  253. func<vector<MatchedTrack>(vector<Seed>,
  254. vector<PixRecHit>,
  255. vector<Track>,
  256. vector<SimHit>,
  257. vector<SimTrack>)>("find_matched_tracks",
  258. FUNC(([](const vector<Seed>& seeds,
  259. const vector<PixRecHit> pixrec_hits,
  260. const vector<Track>& tracks,
  261. const vector<SimHit> sim_hits,
  262. const vector<SimTrack> sim_tracks){
  263. INFO("New event");
  264. INFO(boost::format("Number of tracks is %d, number of seeds is %d, number of hits %d, number of simhits %d")
  265. % tracks.size() % seeds.size() % pixrec_hits.size() % sim_hits.size());
  266. vector<MatchedTrack> matched_tracks;
  267. for(const Track &track : tracks){
  268. const Seed& seed = seeds[track.seedIdx];
  269. if( seed.algoOriginal < 0 || seed.algoOriginal >= seedTypes.size()) continue;
  270. INFO(boost::format(" seed indx= %d algorithm= %d type %s\n")
  271. % track.seedIdx % seed.algoOriginal % seedTypes[seed.algoOriginal]);
  272. for(auto tup : boost::combine(seed.hitIdx, seed.hitType)){ //for hits in track's seed
  273. int hitIdx, hitType;
  274. boost::tie(hitIdx, hitType) = tup;
  275. if(hitType != HIT_TYPE_PIXEL) continue; // take only pixel hits for now
  276. const PixRecHit &rec_hit = pixrec_hits[hitIdx];
  277. TVector3 recoHitPoint;
  278. recoHitPoint.SetXYZ( rec_hit.x, rec_hit.y, rec_hit.z);
  279. float Rreco = recoHitPoint.Perp();
  280. INFO(boost::format(" hit= %d is in %s-L%d at radius= %f\n")
  281. % hitIdx % subdet_names[rec_hit.det] % rec_hit.lay % Rreco);
  282. for(const int& simHitIdx : rec_hit.simHitIdx){ // for sim-hits matching rec-hit
  283. const SimHit &sim_hit = sim_hits[simHitIdx];
  284. const SimTrack &sim_track = sim_tracks[sim_hit.simTrkIdx];
  285. for(const int &trkIdx : sim_track.trkIdx){ // for tracks matched to this sim_track
  286. if(trkIdx == track.idx){ // sim_track matches with our original track
  287. matched_tracks.push_back({seed, rec_hit, track, sim_hit, sim_track});
  288. break;
  289. }
  290. }
  291. }
  292. }
  293. }
  294. return matched_tracks;
  295. })));
  296. fv::apply(find_matched_tracks, fv::tuple(seeds, pixrec_hits, tracks, sim_hits, sim_tracks), "matched_tracks");
  297. }
  298. void setup_residuals(TrackingDataSet &tds){
  299. auto matched_tracks = lookup<vector<MatchedTrack>>("matched_tracks");
  300. // Break matched_tracks into catagories based on Barrel/Endcap and layer
  301. auto matched_tracks_B1 = filter(func<bool(MatchedTrack)>("matched_tracks_B1",
  302. FUNC(([](const MatchedTrack& matched_track){
  303. return in_det(matched_track, PIXEL_BARREL, 1);
  304. }))), matched_tracks);
  305. auto matched_tracks_B2 = filter(func<bool(MatchedTrack)>("matched_tracks_B2",
  306. FUNC(([](const MatchedTrack& matched_track){
  307. return in_det(matched_track, PIXEL_BARREL, 2);
  308. }))), matched_tracks);
  309. auto matched_tracks_B3 = filter(func<bool(MatchedTrack)>("matched_tracks_B3",
  310. FUNC(([](const MatchedTrack& matched_track){
  311. return in_det(matched_track, PIXEL_BARREL, 3);
  312. }))), matched_tracks);
  313. auto matched_tracks_B4 = filter(func<bool(MatchedTrack)>("matched_tracks_B4",
  314. FUNC(([](const MatchedTrack& matched_track){
  315. return in_det(matched_track, PIXEL_BARREL, 4);
  316. }))), matched_tracks);
  317. auto matched_tracks_F1 = filter(func<bool(MatchedTrack)>("matched_tracks_F1",
  318. FUNC(([](const MatchedTrack& matched_track){
  319. return in_det(matched_track, PIXEL_ENDCAP, 1);
  320. }))), matched_tracks);
  321. auto matched_tracks_F2 = filter(func<bool(MatchedTrack)>("matched_tracks_F2",
  322. FUNC(([](const MatchedTrack& matched_track){
  323. return in_det(matched_track, PIXEL_ENDCAP, 2);
  324. }))), matched_tracks);
  325. auto matched_tracks_F3 = filter(func<bool(MatchedTrack)>("matched_tracks_F3",
  326. FUNC(([](const MatchedTrack& matched_track){
  327. return in_det(matched_track, PIXEL_ENDCAP, 3);
  328. }))), matched_tracks);
  329. TH2Params params = {"$\\eta$", 100, -4, 4,
  330. "Residuals(cm)", 100, 0, .01};
  331. auto& calc_residuals_v_eta = func<pair_vec(vector<MatchedTrack>)>("matched_track_hit_residuals",
  332. FUNC(([](const vector<MatchedTrack>& matched_tracks){
  333. vector<float> residuals;
  334. vector<float> etas;
  335. for(auto matched_track : matched_tracks){
  336. auto& pixrec_hit = std::get<PixRecHit>(matched_track);
  337. auto& sim_hit = std::get<SimHit>(matched_track);
  338. auto& sim_track = std::get<SimTrack>(matched_track);
  339. residuals.push_back(displacement(pixrec_hit, sim_hit));
  340. etas.push_back(sim_track.eta);
  341. }
  342. return std::make_pair(etas, residuals);
  343. })));
  344. tds.register_container<ContainerTH2Many<float>>("matched_track_residuals_v_eta_all",
  345. fv::apply(calc_residuals_v_eta, fv::tuple(matched_tracks)),
  346. "Matched Track Residuals - All", params);
  347. tds.register_container<ContainerTH2Many<float>>("matched_track_residuals_v_eta_B1",
  348. fv::apply(calc_residuals_v_eta, fv::tuple(matched_tracks_B1)),
  349. "Matched Track Residuals - B1", params);
  350. tds.register_container<ContainerTH2Many<float>>("matched_track_residuals_v_eta_B2",
  351. fv::apply(calc_residuals_v_eta, fv::tuple(matched_tracks_B2)),
  352. "Matched Track Residuals - B2", params);
  353. tds.register_container<ContainerTH2Many<float>>("matched_track_residuals_v_eta_B3",
  354. fv::apply(calc_residuals_v_eta, fv::tuple(matched_tracks_B3)),
  355. "Matched Track Residuals - B3", params);
  356. tds.register_container<ContainerTH2Many<float>>("matched_track_residuals_v_eta_B4",
  357. fv::apply(calc_residuals_v_eta, fv::tuple(matched_tracks_B4)),
  358. "Matched Track Residuals - B4", params);
  359. tds.register_container<ContainerTH2Many<float>>("matched_track_residuals_v_eta_F1",
  360. fv::apply(calc_residuals_v_eta, fv::tuple(matched_tracks_F1)),
  361. "Matched Track Residuals - F1", params);
  362. tds.register_container<ContainerTH2Many<float>>("matched_track_residuals_v_eta_F2",
  363. fv::apply(calc_residuals_v_eta, fv::tuple(matched_tracks_F2)),
  364. "Matched Track Residuals - F2", params);
  365. tds.register_container<ContainerTH2Many<float>>("matched_track_residuals_v_eta_F3",
  366. fv::apply(calc_residuals_v_eta, fv::tuple(matched_tracks_F3)),
  367. "Matched Track Residuals - F3", params);
  368. params = {"$\\eta$", 100, -4, 4,
  369. "$\\Delta \\phi$(rad)", 100, -.002, .002};
  370. auto& calc_dphi_v_eta = func<pair_vec(vector<MatchedTrack>)>("matched_track_hit_dphis",
  371. FUNC(([](const vector<MatchedTrack>& matched_tracks){
  372. vector<float> dphis;
  373. vector<float> etas;
  374. for(auto matched_track : matched_tracks){
  375. auto& pixrec_hit = std::get<PixRecHit>(matched_track);
  376. auto& sim_hit = std::get<SimHit>(matched_track);
  377. auto& sim_track = std::get<SimTrack>(matched_track);
  378. dphis.push_back(atan2(sim_hit.x, sim_hit.y) - atan2(pixrec_hit.x, pixrec_hit.y));
  379. etas.push_back(sim_track.eta);
  380. }
  381. return std::make_pair(etas, dphis);
  382. })));
  383. tds.register_container<ContainerTH2Many<float>>("matched_track_dphis_v_eta_all",
  384. fv::apply(calc_dphi_v_eta, fv::tuple(matched_tracks)),
  385. "Matched Track $\\Delta \\phi$ - All", params);
  386. tds.register_container<ContainerTH2Many<float>>("matched_track_dphis_v_eta_B1",
  387. fv::apply(calc_dphi_v_eta, fv::tuple(matched_tracks_B1)),
  388. "Matched Track $\\Delta \\phi$ - B1", params);
  389. tds.register_container<ContainerTH2Many<float>>("matched_track_dphis_v_eta_B2",
  390. fv::apply(calc_dphi_v_eta, fv::tuple(matched_tracks_B2)),
  391. "Matched Track $\\Delta \\phi$ - B2", params);
  392. tds.register_container<ContainerTH2Many<float>>("matched_track_dphis_v_eta_B3",
  393. fv::apply(calc_dphi_v_eta, fv::tuple(matched_tracks_B3)),
  394. "Matched Track $\\Delta \\phi$ - B3", params);
  395. tds.register_container<ContainerTH2Many<float>>("matched_track_dphis_v_eta_B4",
  396. fv::apply(calc_dphi_v_eta, fv::tuple(matched_tracks_B4)),
  397. "Matched Track $\\Delta \\phi$ - B4", params);
  398. tds.register_container<ContainerTH2Many<float>>("matched_track_dphis_v_eta_F1",
  399. fv::apply(calc_dphi_v_eta, fv::tuple(matched_tracks_F1)),
  400. "Matched Track $\\Delta \\phi$ - F1", params);
  401. tds.register_container<ContainerTH2Many<float>>("matched_track_dphis_v_eta_F2",
  402. fv::apply(calc_dphi_v_eta, fv::tuple(matched_tracks_F2)),
  403. "Matched Track $\\Delta \\phi$ - F2", params);
  404. tds.register_container<ContainerTH2Many<float>>("matched_track_dphis_v_eta_F3",
  405. fv::apply(calc_dphi_v_eta, fv::tuple(matched_tracks_F3)),
  406. "Matched Track $\\Delta \\phi$ - F3", params);
  407. auto& calc_dz_v_eta = func<pair_vec(vector<MatchedTrack>)>("matched_track_hit_dphis",
  408. FUNC(([](const vector<MatchedTrack>& matched_tracks){
  409. vector<float> dzs;
  410. vector<float> etas;
  411. for(auto matched_track : matched_tracks){
  412. auto& pixrec_hit = std::get<PixRecHit>(matched_track);
  413. auto& sim_hit = std::get<SimHit>(matched_track);
  414. auto& sim_track = std::get<SimTrack>(matched_track);
  415. dzs.push_back(sim_hit.z - pixrec_hit.z);
  416. etas.push_back(sim_track.eta);
  417. }
  418. return std::make_pair(etas, dzs);
  419. })));
  420. params = {"$\\eta$", 100, -4, 4,
  421. "$\\Delta z$(cm)", 100, -.002, .002};
  422. tds.register_container<ContainerTH2Many<float>>("matched_track_dzs_v_eta_all",
  423. fv::apply(calc_dz_v_eta, fv::tuple(matched_tracks)),
  424. "Matched Track $\\Delta z$ - All", params);
  425. tds.register_container<ContainerTH2Many<float>>("matched_track_dzs_v_eta_B1",
  426. fv::apply(calc_dz_v_eta, fv::tuple(matched_tracks_B1)),
  427. "Matched Track $\\Delta z$ - B1", params);
  428. tds.register_container<ContainerTH2Many<float>>("matched_track_dzs_v_eta_B2",
  429. fv::apply(calc_dz_v_eta, fv::tuple(matched_tracks_B2)),
  430. "Matched Track $\\Delta z$ - B2", params);
  431. tds.register_container<ContainerTH2Many<float>>("matched_track_dzs_v_eta_B3",
  432. fv::apply(calc_dz_v_eta, fv::tuple(matched_tracks_B3)),
  433. "Matched Track $\\Delta z$ - B3", params);
  434. tds.register_container<ContainerTH2Many<float>>("matched_track_dzs_v_eta_B4",
  435. fv::apply(calc_dz_v_eta, fv::tuple(matched_tracks_B4)),
  436. "Matched Track $\\Delta z$ - B4", params);
  437. auto& calc_dr_v_eta = func<pair_vec(vector<MatchedTrack>)>("matched_track_hit_drs",
  438. FUNC(([](const vector<MatchedTrack>& matched_tracks){
  439. vector<float> drs;
  440. vector<float> etas;
  441. for(auto matched_track : matched_tracks){
  442. auto& rec_hit = std::get<PixRecHit>(matched_track);
  443. auto& sim_hit = std::get<SimHit>(matched_track);
  444. float sim_hit_r = sqrt(pow(sim_hit.x, 2) + pow(sim_hit.y, 2));
  445. float rec_hit_r = sqrt(pow(rec_hit.x, 2) + pow(rec_hit.y, 2));
  446. drs.push_back(sim_hit_r - rec_hit_r);
  447. etas.push_back(std::get<SimTrack>(matched_track).eta);
  448. }
  449. return std::make_pair(etas, drs);
  450. })));
  451. params = {"$\\eta$", 100, -4, 4,
  452. "$\\Delta r$(cm)", 100, -.002, .002};
  453. tds.register_container<ContainerTH2Many<float>>("matched_track_drs_v_eta_F1",
  454. fv::apply(calc_dr_v_eta, fv::tuple(matched_tracks_F1)),
  455. "Matched Track $\\Delta r$ - F1", params);
  456. tds.register_container<ContainerTH2Many<float>>("matched_track_drs_v_eta_F2",
  457. fv::apply(calc_dr_v_eta, fv::tuple(matched_tracks_F2)),
  458. "Matched Track $\\Delta r$ - F2", params);
  459. tds.register_container<ContainerTH2Many<float>>("matched_track_drs_v_eta_F3",
  460. fv::apply(calc_dr_v_eta, fv::tuple(matched_tracks_F3)),
  461. "Matched Track $\\Delta r$ - F3", params);
  462. }
  463. void setup_ladder_vs_phi(TrackingDataSet& tds){
  464. auto& calc_phi_vs_ladder =
  465. func<pair_vec(vector<PixRecHit>, int)>("calc_phi_vs_ladder",
  466. FUNC(([](const vector<PixRecHit>& pixrec_hits, int layer){
  467. std::vector<float> phis;
  468. std::vector<float> ladders;
  469. for(const PixRecHit &pixrec_hit : pixrec_hits){
  470. if(pixrec_hit.det == PIXEL_BARREL && pixrec_hit.lay == layer){
  471. float phi = pixrec_hit.phi();
  472. unsigned short ladder = pixrec_hit.ladder_blade;
  473. phis.push_back(phi);
  474. ladders.push_back(ladder);
  475. INFO(boost::format("Hit Phi: %7.3f , Hit detID: %8x, Ladder: %d")
  476. % phi % pixrec_hit.detId % ladder);
  477. }
  478. }
  479. return std::make_pair(phis, ladders);
  480. })));
  481. TH2Params params = {"$\\phi$", 100, -3.14159, 3.14159,
  482. "Ladder Number", 50, 0, 50};
  483. tds.register_container<ContainerTH2Many<float>>("rechit_phi_vs_ladder_l1",
  484. fv::apply(calc_phi_vs_ladder, fv::tuple(pixrec_hits, constant<int>("1",1))),
  485. "phi vs ladder number - Layer 1", params);
  486. tds.register_container<ContainerTH2Many<float>>("rechit_phi_vs_ladder_l2",
  487. fv::apply(calc_phi_vs_ladder, fv::tuple(pixrec_hits, constant<int>("2",2))),
  488. "phi vs ladder number - Layer 2", params);
  489. tds.register_container<ContainerTH2Many<float>>("rechit_phi_vs_ladder_l3",
  490. fv::apply(calc_phi_vs_ladder, fv::tuple(pixrec_hits, constant<int>("3",3))),
  491. "phi vs ladder number - Layer 3", params);
  492. }
  493. void run_analysis(const vector<fv::util::DataFileDescriptor>& dfds, const string output_filename, bool silent){
  494. gSystem->Load("libfilval.so");
  495. auto replace_suffix = [](const std::string& input, const std::string& new_suffix){
  496. return input.substr(0, input.find_last_of(".")) + new_suffix;
  497. };
  498. string log_filename = replace_suffix(output_filename, ".log");
  499. fv::util::Log::init_logger(log_filename, fv::util::LogPriority::kLogInfo);
  500. TrackingDataSet tds(output_filename, dfds, "trackingNtuple/tree");
  501. /* tds.set_max_events(10); */
  502. register_objects(tds);
  503. /* setup_matched_tracks(tds); */
  504. /* setup_residuals(tds); */
  505. setup_first_hit_pairs(tds);
  506. /* setup_ladder_vs_phi(tds); */
  507. tds.process(silent);
  508. tds.save_all();
  509. }
  510. int main(int argc, char * argv[]){
  511. fv::util::ArgParser args(argc, argv);
  512. bool silent = args.cmdOptionExists("-s");
  513. string output_filename = args.cmdOptionExists("-o") ? args.getCmdOption("-o") : "output.root";
  514. if(args.cmdOptionExists("-F")){
  515. auto file_list = fv::util::read_input_list(args.getCmdOption("-F"));
  516. run_analysis(file_list, output_filename, silent);
  517. }else if(args.cmdOptionExists("-f")){
  518. string input_filename = args.getCmdOption("-f");
  519. string data_label = args.cmdOptionExists("-l") ? args.getCmdOption("-l") : "";
  520. fv::util::DataFileDescriptor dfd(input_filename, data_label);
  521. run_analysis(std::vector<fv::util::DataFileDescriptor>({dfd}), output_filename, silent);
  522. }else {
  523. cout << "Usage: ./tracking_validation (-s) {-o output_filename} -F datafiles.txt" << endl;
  524. cout << " ./tracking_validation (-s) {-l DATA_LABEL} {-o output_filename} -f treefile.root" << endl;
  525. }
  526. }