tracking_validation.cpp 37 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647
  1. #include <iostream>
  2. #include <vector>
  3. #include <map>
  4. #include <utility>
  5. #include <numeric>
  6. #include <limits>
  7. #include <TSystem.h>
  8. #include "filval/filval.hpp"
  9. #include "filval/root/filval.hpp"
  10. #include <boost/range/combine.hpp>
  11. #include <boost/format.hpp>
  12. #include "analysis/TrackingNtuple.h"
  13. #include "analysis/obj_types.cpp"
  14. #include "analysis/common.hpp"
  15. using namespace std;
  16. using namespace fv;
  17. using namespace fv::root;
  18. typedef std::pair<std::vector<float>,std::vector<float>> pair_vec;
  19. #define HIT_TYPE_PIXEL 0
  20. #define HIT_TYPE_GLUED 1
  21. #define HIT_TYPE_STRIP 2
  22. #define PIXEL_BARREL 1
  23. #define PIXEL_ENDCAP 2
  24. std::map<int,string> subdet_names = {{1, "BPIX"}, {2, "FPIX"}};
  25. typedef std::tuple<PixRecHit, SimHit> HitPair;
  26. void setup_hit_pair_functions(){
  27. func<pair_vec(vector<HitPair>)>("calc_dphi_v_eta",
  28. FUNC(([](const vector<HitPair>& hit_pairs){
  29. vector<float> dphis;
  30. vector<float> etas;
  31. for(auto hit_pair : hit_pairs){
  32. auto& pixrec_hit = std::get<PixRecHit>(hit_pair);
  33. auto& sim_hit = std::get<SimHit>(hit_pair);
  34. dphis.push_back(atan2(sim_hit.x, sim_hit.y) - atan2(pixrec_hit.x, pixrec_hit.y));
  35. etas.push_back(pseudorapidity(pixrec_hit));
  36. }
  37. return std::make_pair(etas, dphis);
  38. })));
  39. func<pair_vec(vector<HitPair>)>("calc_dr_v_eta",
  40. FUNC(([](const vector<HitPair>& hit_pairs){
  41. vector<float> drs;
  42. vector<float> etas;
  43. for(auto hit_pair : hit_pairs){
  44. auto& pixrec_hit = std::get<PixRecHit>(hit_pair);
  45. auto& sim_hit = std::get<SimHit>(hit_pair);
  46. drs.push_back(displacement(pixrec_hit, sim_hit));
  47. etas.push_back(pseudorapidity(pixrec_hit));
  48. }
  49. return std::make_pair(etas, drs);
  50. })));
  51. func<pair_vec(vector<HitPair>)>("calc_dz_v_eta",
  52. FUNC(([](const vector<HitPair>& hit_pairs){
  53. vector<float> dzs;
  54. vector<float> etas;
  55. for(auto hit_pair : hit_pairs){
  56. auto& pixrec_hit = std::get<PixRecHit>(hit_pair);
  57. auto& sim_hit = std::get<SimHit>(hit_pair);
  58. dzs.push_back(sim_hit.z - pixrec_hit.z);
  59. etas.push_back(pseudorapidity(pixrec_hit));
  60. }
  61. return std::make_pair(etas, dzs);
  62. })));
  63. func<pair_vec(vector<HitPair>)>("calc_drho_v_eta",
  64. FUNC(([](const vector<HitPair>& hit_pairs){
  65. vector<float> drhos;
  66. vector<float> etas;
  67. for(auto hit_pair : hit_pairs){
  68. auto& pixrec_hit = std::get<PixRecHit>(hit_pair);
  69. auto& sim_hit = std::get<SimHit>(hit_pair);
  70. drhos.push_back(rho(sim_hit) - rho(pixrec_hit));
  71. etas.push_back(pseudorapidity(pixrec_hit));
  72. }
  73. return std::make_pair(etas, drhos);
  74. })));
  75. }
  76. vector<string> seedTypes =
  77. {"initialStepSeeds",
  78. "highPtTripletStepSeeds",
  79. "mixedTripletStepSeeds",
  80. "pixelLessStepSeeds",
  81. "tripletElectronSeeds",
  82. "pixelPairElectronSeeds",
  83. "stripPairElectronSeeds"};
  84. Value<vector<SuperCluster>>* super_clusters;
  85. Value<vector<Seed>>* seeds;
  86. Value<vector<PixRecHit>>* pixrec_hits;
  87. Value<vector<Track>>* tracks;
  88. Value<vector<SimHit>>* sim_hits;
  89. Value<vector<SimTrack>>* sim_tracks;
  90. void register_objects(TrackingDataSet& tds){
  91. super_clusters = register_super_clusters(tds);
  92. seeds = register_seeds(tds);
  93. pixrec_hits = register_pixrec_hits(tds);
  94. tracks = register_tracks(tds);
  95. sim_hits = register_sim_hits(tds);
  96. sim_tracks = register_sim_tracks(tds);
  97. }
  98. template<typename A>
  99. bool in_det(const A& hit, const unsigned int& det, const unsigned int& lay){
  100. return hit.det == det && hit.lay == lay;
  101. }
  102. void setup_skipped_layer_hit_pairs(TrackingDataSet& tds){
  103. /* Finds SimHit/RecHit pairs
  104. */
  105. auto find_matched_nth_hit_in_layer_with_skip =
  106. func<vector<HitPair>(vector<Seed>,
  107. vector<PixRecHit>,
  108. vector<SimHit>,
  109. vector<Track>,
  110. vector<SimTrack>,
  111. int, int, int)>("find_matched_nth_hit_in_layer_with_skip",
  112. FUNC(([](const vector<Seed>& seeds,
  113. const vector<PixRecHit>& pixrec_hits,
  114. const vector<SimHit>& sim_hits,
  115. const vector<Track>& tracks,
  116. const vector<SimTrack>& sim_tracks,
  117. const unsigned int&& det,
  118. const unsigned int&& pix_layer,
  119. const unsigned int&& skip){
  120. vector<HitPair> matched_hits;
  121. for(const Track &trk : tracks){ // loop over all tracks
  122. const Seed &seed = seeds[trk.seedIdx];
  123. // Require seed's source algorithm is one of those in seedTypes
  124. if(seed.algoOriginal < 0 || seed.algoOriginal >= seedTypes.size()) continue;
  125. // Require seed w/ at least 2 hits
  126. if(seed.hitIdx.size() < 2) continue;
  127. // take only pixel hits for now
  128. if(seed.hitType[0] != HIT_TYPE_PIXEL || seed.hitType[1] != HIT_TYPE_PIXEL) continue;
  129. const PixRecHit &rec_hit1 = pixrec_hits[seed.hitIdx[0]];
  130. const PixRecHit &rec_hit2 = pixrec_hits[seed.hitIdx[1]];
  131. if(in_det(rec_hit1, det, pix_layer) && in_det(rec_hit2, det, pix_layer+1+skip)){
  132. // We have the RecHit we want to work with, now find a properly* matched SimHit
  133. if(rec_hit1.simHitIdx.size() == 0) continue; // if rechit is matched to no simhits, give up.
  134. if(trk.simTrkIdx.size() == 0) continue; // if track is matched to no simtracks, give up.
  135. const int &simTrkIdx = trk.simTrkIdx[0]; // get first matched track
  136. for(const int& simHitIdx : sim_tracks[simTrkIdx].simHitIdx){ // loop over SimHits from SimTrack
  137. if(simHitIdx == rec_hit1.simHitIdx[0]) // take first matched simhit (should be the closest one)
  138. matched_hits.push_back({rec_hit1, sim_hits[rec_hit1.simHitIdx[0]]});
  139. }
  140. }
  141. }
  142. return matched_hits;
  143. })));
  144. auto pix_barrel = constant("PIXEL_BARREL", PIXEL_BARREL);
  145. auto pix_endcap = constant("PIXEL_ENDCAP", PIXEL_ENDCAP);
  146. auto L1 = constant("L1", 1);
  147. auto L2 = constant("L2", 2);
  148. auto skip_zero = constant("skip_zero", 0);
  149. auto skip_one = constant("skip_one", 1);
  150. // First hit on 1st/2nd bpix layers and second hit in 2nd/3rd
  151. auto first_hits_in_B1_skip_zero = fv::tup_apply(find_matched_nth_hit_in_layer_with_skip,
  152. fv::tuple(seeds, pixrec_hits, sim_hits, tracks, sim_tracks, pix_barrel, L1, skip_zero), "first_hits_in_B1_skip_zero");
  153. auto first_hits_in_B2_skip_zero = fv::tup_apply(find_matched_nth_hit_in_layer_with_skip,
  154. fv::tuple(seeds, pixrec_hits, sim_hits, tracks, sim_tracks, pix_barrel, L2, skip_zero), "first_hits_in_B2_skip_zero");
  155. // First hit on 1st/2nd fpix layers and second hit in 2nd/3rd
  156. auto first_hits_in_F1_skip_zero = fv::tup_apply(find_matched_nth_hit_in_layer_with_skip,
  157. fv::tuple(seeds, pixrec_hits, sim_hits, tracks, sim_tracks, pix_endcap, L1, skip_zero), "first_hits_in_F1_skip_zero");
  158. auto first_hits_in_F2_skip_zero = fv::tup_apply(find_matched_nth_hit_in_layer_with_skip,
  159. fv::tuple(seeds, pixrec_hits, sim_hits, tracks, sim_tracks, pix_endcap, L2, skip_zero), "first_hits_in_F2_skip_zero");
  160. // First hit on 1st/2nd fpix layers and second hit in 3rd/4th
  161. auto first_hits_in_B1_skip_one = fv::tup_apply(find_matched_nth_hit_in_layer_with_skip,
  162. fv::tuple(seeds, pixrec_hits, sim_hits, tracks, sim_tracks, pix_barrel, L1, skip_one), "first_hits_in_B1_skip_one");
  163. auto first_hits_in_B2_skip_one = fv::tup_apply(find_matched_nth_hit_in_layer_with_skip,
  164. fv::tuple(seeds, pixrec_hits, sim_hits, tracks, sim_tracks, pix_barrel, L2, skip_one), "first_hits_in_B2_skip_one");
  165. // First hit on 1st fpix layer and second hit in 3rd layer
  166. auto first_hits_in_F1_skip_one = fv::tup_apply(find_matched_nth_hit_in_layer_with_skip,
  167. fv::tuple(seeds, pixrec_hits, sim_hits, tracks, sim_tracks, pix_endcap, L1, skip_one), "first_hits_in_F1_skip_one");
  168. TH2Params params_dz = {"$\\eta$", 100, -4, 4,
  169. "$\\Delta z$(rad)", 50, -.01, .01};
  170. auto calc_dz_v_eta = lookup_function<pair_vec(vector<HitPair>)>("calc_dz_v_eta");
  171. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_first_hits_in_B1_skip_zero",
  172. fv::apply(calc_dz_v_eta, first_hits_in_B1_skip_zero),
  173. "First Hit in BPIX-L1 - Skip 0", params_dz);
  174. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_first_hits_in_B2_skip_zero",
  175. fv::apply(calc_dz_v_eta, first_hits_in_B2_skip_zero),
  176. "First Hit in BPIX-L2 - Skip 0", params_dz);
  177. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_first_hits_in_F1_skip_zero",
  178. fv::apply(calc_dz_v_eta, first_hits_in_F1_skip_zero),
  179. "First Hit in FPIX-L1 - Skip 0", params_dz);
  180. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_first_hits_in_F2_skip_zero",
  181. fv::apply(calc_dz_v_eta, first_hits_in_F2_skip_zero),
  182. "First Hit in FPIX-L2 - Skip 0", params_dz);
  183. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_first_hits_in_B1_skip_one",
  184. fv::apply(calc_dz_v_eta, first_hits_in_B1_skip_one),
  185. "First Hit in BPIX-L1 - Skip 1", params_dz);
  186. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_first_hits_in_B2_skip_one",
  187. fv::apply(calc_dz_v_eta, first_hits_in_B2_skip_one),
  188. "First Hit in BPIX-L2 - Skip 1", params_dz);
  189. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_first_hits_in_F1_skip_one",
  190. fv::apply(calc_dz_v_eta, first_hits_in_F1_skip_one),
  191. "First Hit in FPIX-L1 - Skip 1", params_dz);
  192. }
  193. void setup_first_hit_pairs(TrackingDataSet& tds){
  194. auto find_matched_nth_hit_in_layer =
  195. func<vector<HitPair>(vector<Seed>,
  196. vector<PixRecHit>,
  197. vector<SimHit>,
  198. vector<Track>,
  199. vector<SimTrack>,
  200. int, int, int)>("find_matched_nth_hit_in_layer",
  201. FUNC(([](const vector<Seed>& seeds,
  202. const vector<PixRecHit>& pixrec_hits,
  203. const vector<SimHit>& sim_hits,
  204. const vector<Track>& tracks,
  205. const vector<SimTrack>& sim_tracks,
  206. const unsigned int&& det,
  207. const unsigned int&& pix_layer,
  208. const unsigned int&& hit_number){
  209. vector<HitPair> matched_hits;
  210. for(const Track &trk : tracks){ // loop over all tracks
  211. const Seed &seed = seeds[trk.seedIdx];
  212. if(seed.hitIdx.size() <= hit_number) continue; // looking for hit_number'th hit
  213. if(seed.algoOriginal < 0 || seed.algoOriginal >= seedTypes.size()) continue;
  214. if(seed.hitType[hit_number] != HIT_TYPE_PIXEL) continue; // take only pixel hits for now
  215. const PixRecHit &rec_hit = pixrec_hits[seed.hitIdx[hit_number]];
  216. if(rec_hit.det == det && rec_hit.lay == pix_layer){
  217. // We have the RecHit we want to work with, now find a properly* matched SimHit
  218. if(rec_hit.simHitIdx.size() == 0) continue; // if rechit is matched to no simhits, give up.
  219. if(trk.simTrkIdx.size() == 0) continue; // if rechit is matched to no simhits, give up.
  220. const int &simTrkIdx = trk.simTrkIdx[0]; // get first matched track
  221. for(const int& simHitIdx : sim_tracks[simTrkIdx].simHitIdx){ // loop over SimHits from SimTrack
  222. if(simHitIdx == rec_hit.simHitIdx[0]) // take first matched simhit (should be the closest one)
  223. matched_hits.push_back({rec_hit, sim_hits[rec_hit.simHitIdx[0]]});
  224. }
  225. }
  226. }
  227. return matched_hits;
  228. })));
  229. auto barrel_val = constant("PIXEL_BARREL", PIXEL_BARREL);
  230. auto endcap_val = constant("PIXEL_ENDCAP", PIXEL_ENDCAP);
  231. auto first_hit = constant("1st", 0);
  232. auto second_hit = constant("2nd", 1);
  233. // First hits on inner three bpix layers
  234. auto first_hits_in_B1 = fv::tup_apply(find_matched_nth_hit_in_layer,
  235. fv::tuple(seeds, pixrec_hits, sim_hits, tracks, sim_tracks, barrel_val, constant("L1", 1), first_hit), "first_hits_in_B1");
  236. auto first_hits_in_B2 = fv::tup_apply(find_matched_nth_hit_in_layer,
  237. fv::tuple(seeds, pixrec_hits, sim_hits, tracks, sim_tracks, barrel_val, constant("L2", 2), first_hit), "first_hits_in_B2");
  238. // Second hits on outer three bpix layers
  239. auto second_hits_in_B2 = fv::tup_apply(find_matched_nth_hit_in_layer,
  240. fv::tuple(seeds, pixrec_hits, sim_hits, tracks, sim_tracks, barrel_val, constant("L2", 2), second_hit), "second_hits_in_B2");
  241. auto second_hits_in_B3 = fv::tup_apply(find_matched_nth_hit_in_layer,
  242. fv::tuple(seeds, pixrec_hits, sim_hits, tracks, sim_tracks, barrel_val, constant("L3", 3), second_hit), "second_hits_in_B3");
  243. auto second_hits_in_B4 = fv::tup_apply(find_matched_nth_hit_in_layer,
  244. fv::tuple(seeds, pixrec_hits, sim_hits, tracks, sim_tracks, barrel_val, constant("L4", 4), second_hit), "second_hits_in_B4");
  245. auto first_hits_in_F1 = fv::tup_apply(find_matched_nth_hit_in_layer,
  246. fv::tuple(seeds, pixrec_hits, sim_hits, tracks, sim_tracks, endcap_val, constant("L1", 1), first_hit), "first_hits_in_F1");
  247. auto first_hits_in_F2 = fv::tup_apply(find_matched_nth_hit_in_layer,
  248. fv::tuple(seeds, pixrec_hits, sim_hits, tracks, sim_tracks, endcap_val, constant("L2", 2), first_hit), "first_hits_in_F2");
  249. auto second_hits_in_F2 = fv::tup_apply(find_matched_nth_hit_in_layer,
  250. fv::tuple(seeds, pixrec_hits, sim_hits, tracks, sim_tracks, endcap_val, constant("L2", 2), second_hit), "second_hits_in_F2");
  251. auto second_hits_in_F3 = fv::tup_apply(find_matched_nth_hit_in_layer,
  252. fv::tuple(seeds, pixrec_hits, sim_hits, tracks, sim_tracks, endcap_val, constant("L3", 3), second_hit), "second_hits_in_F3");
  253. // Even vs Odd Ladders
  254. auto even = constant("even", false);
  255. auto odd = constant("odd", true);
  256. auto select_even_odd_ladder_blade_hit_pairs =
  257. func<vector<HitPair>(vector<HitPair>, bool)>("select_even_odd_ladder_blade_hit_pairs",
  258. FUNC(([](const vector<HitPair>& hit_pairs,
  259. const bool &&odd){
  260. vector<HitPair> even_pairs;
  261. for(const HitPair &hit_pair : hit_pairs){ // loop over all seeds
  262. if(std::get<PixRecHit>(hit_pair).ladder_blade % 2 == odd){
  263. even_pairs.push_back(hit_pair);
  264. }
  265. }
  266. return even_pairs;
  267. })));
  268. auto first_hits_in_B1_even_ladder = fv::tup_apply(select_even_odd_ladder_blade_hit_pairs,
  269. fv::tuple(first_hits_in_B1, even), "first_hits_in_B1_even_ladder");
  270. auto first_hits_in_B1_odd_ladder = fv::tup_apply(select_even_odd_ladder_blade_hit_pairs,
  271. fv::tuple(first_hits_in_B1, odd), "first_hits_in_B1_odd_ladder");
  272. auto first_hits_in_B2_even_ladder = fv::tup_apply(select_even_odd_ladder_blade_hit_pairs,
  273. fv::tuple(first_hits_in_B2, even), "first_hits_in_B2_even_ladder");
  274. auto first_hits_in_B2_odd_ladder = fv::tup_apply(select_even_odd_ladder_blade_hit_pairs,
  275. fv::tuple(first_hits_in_B2, odd), "first_hits_in_B2_odd_ladder");
  276. //Plots for dPhi of collections defined above
  277. auto calc_dphi_v_eta = lookup_function<pair_vec(vector<HitPair>)>("calc_dphi_v_eta");
  278. TH2Params params_dphi = {"$\\eta$", 100, -4, 4,
  279. "$\\Delta \\phi$(rad)", 50, -.0015, .0015};
  280. tds.register_container<ContainerTH2Many<float>>("dphi_v_eta_first_hits_in_B1",
  281. fv::apply(calc_dphi_v_eta, first_hits_in_B1),
  282. "First Hit in BPIX-L1", params_dphi);
  283. tds.register_container<ContainerTH2Many<float>>("dphi_v_eta_first_hits_in_B2",
  284. fv::apply(calc_dphi_v_eta, first_hits_in_B2),
  285. "First Hit in BPIX-L2", params_dphi);
  286. tds.register_container<ContainerTH2Many<float>>("dphi_v_eta_first_hits_in_B1_even_ladder",
  287. fv::apply(calc_dphi_v_eta, first_hits_in_B1_even_ladder),
  288. "First Hit in BPIX-L1 - Even Ladders", params_dphi);
  289. tds.register_container<ContainerTH2Many<float>>("dphi_v_eta_first_hits_in_B1_odd_ladder",
  290. fv::apply(calc_dphi_v_eta, first_hits_in_B1_odd_ladder),
  291. "First Hit in BPIX-L1 - Odd Ladders", params_dphi);
  292. tds.register_container<ContainerTH2Many<float>>("dphi_v_eta_first_hits_in_B2_even_ladder",
  293. fv::apply(calc_dphi_v_eta, first_hits_in_B2_even_ladder),
  294. "First Hit in BPIX-L2 - Even Ladders", params_dphi);
  295. tds.register_container<ContainerTH2Many<float>>("dphi_v_eta_first_hits_in_B2_odd_ladder",
  296. fv::apply(calc_dphi_v_eta, first_hits_in_B2_odd_ladder),
  297. "First Hit in BPIX-L2 - Odd Ladders", params_dphi);
  298. tds.register_container<ContainerTH2Many<float>>("dphi_v_eta_second_hits_in_B2",
  299. fv::apply(calc_dphi_v_eta, second_hits_in_B2),
  300. "Second Hit in BPIX-L2", params_dphi);
  301. tds.register_container<ContainerTH2Many<float>>("dphi_v_eta_second_hits_in_B3",
  302. fv::apply(calc_dphi_v_eta, second_hits_in_B3),
  303. "Second Hit in BPIX-L3", params_dphi);
  304. tds.register_container<ContainerTH2Many<float>>("dphi_v_eta_second_hits_in_B4",
  305. fv::apply(calc_dphi_v_eta, second_hits_in_B4),
  306. "Second Hit in BPIX-L4", params_dphi);
  307. tds.register_container<ContainerTH2Many<float>>("dphi_v_eta_first_hits_in_F1",
  308. fv::apply(calc_dphi_v_eta, first_hits_in_F1),
  309. "First Hit in FPIX-L1", params_dphi);
  310. tds.register_container<ContainerTH2Many<float>>("dphi_v_eta_first_hits_in_F2",
  311. fv::apply(calc_dphi_v_eta, first_hits_in_F2),
  312. "First Hit in FPIX-L2", params_dphi);
  313. tds.register_container<ContainerTH2Many<float>>("dphi_v_eta_second_hits_in_F2",
  314. fv::apply(calc_dphi_v_eta, second_hits_in_F2),
  315. "Second Hit in FPIX-L2", params_dphi);
  316. tds.register_container<ContainerTH2Many<float>>("dphi_v_eta_second_hits_in_F3",
  317. fv::apply(calc_dphi_v_eta, second_hits_in_F3),
  318. "Second Hit in FPIX-L3", params_dphi);
  319. //Plots for dr of collections defined above
  320. auto calc_dr_v_eta = lookup_function<pair_vec(vector<HitPair>)>("calc_dr_v_eta");
  321. TH2Params params_dr = {"$\\eta$", 100, -4, 4,
  322. "$\\Delta \\r$(rad)", 50, 0, .006};
  323. tds.register_container<ContainerTH2Many<float>>("dr_v_eta_first_hits_in_B1",
  324. fv::apply(calc_dr_v_eta, first_hits_in_B1),
  325. "First Hit in BPIX-L1", params_dr);
  326. tds.register_container<ContainerTH2Many<float>>("dr_v_eta_first_hits_in_B2",
  327. fv::apply(calc_dr_v_eta, first_hits_in_B2),
  328. "First Hit in BPIX-L2", params_dr);
  329. tds.register_container<ContainerTH2Many<float>>("dr_v_eta_first_hits_in_B1_even_ladder",
  330. fv::apply(calc_dr_v_eta, first_hits_in_B1_even_ladder),
  331. "First Hit in BPIX-L1 - Even Ladders", params_dr);
  332. tds.register_container<ContainerTH2Many<float>>("dr_v_eta_first_hits_in_B1_odd_ladder",
  333. fv::apply(calc_dr_v_eta, first_hits_in_B1_odd_ladder),
  334. "First Hit in BPIX-L1 - Odd Ladders", params_dr);
  335. tds.register_container<ContainerTH2Many<float>>("dr_v_eta_first_hits_in_B2_even_ladder",
  336. fv::apply(calc_dr_v_eta, first_hits_in_B2_even_ladder),
  337. "First Hit in BPIX-L2 - Even Ladders", params_dr);
  338. tds.register_container<ContainerTH2Many<float>>("dr_v_eta_first_hits_in_B2_odd_ladder",
  339. fv::apply(calc_dr_v_eta, first_hits_in_B2_odd_ladder),
  340. "First Hit in BPIX-L2 - Odd Ladders", params_dr);
  341. tds.register_container<ContainerTH2Many<float>>("dr_v_eta_second_hits_in_B2",
  342. fv::apply(calc_dr_v_eta, second_hits_in_B2),
  343. "Second Hit in BPIX-L2", params_dr);
  344. tds.register_container<ContainerTH2Many<float>>("dr_v_eta_second_hits_in_B3",
  345. fv::apply(calc_dr_v_eta, second_hits_in_B3),
  346. "Second Hit in BPIX-L3", params_dr);
  347. tds.register_container<ContainerTH2Many<float>>("dr_v_eta_second_hits_in_B4",
  348. fv::apply(calc_dr_v_eta, second_hits_in_B4),
  349. "Second Hit in BPIX-L4", params_dr);
  350. tds.register_container<ContainerTH2Many<float>>("dr_v_eta_first_hits_in_F1",
  351. fv::apply(calc_dr_v_eta, first_hits_in_F1),
  352. "First Hit in FPIX-L1", params_dr);
  353. tds.register_container<ContainerTH2Many<float>>("dr_v_eta_first_hits_in_F2",
  354. fv::apply(calc_dr_v_eta, first_hits_in_F2),
  355. "First Hit in FPIX-L2", params_dr);
  356. tds.register_container<ContainerTH2Many<float>>("dr_v_eta_second_hits_in_F2",
  357. fv::apply(calc_dr_v_eta, second_hits_in_F2),
  358. "Second Hit in FPIX-L2", params_dr);
  359. tds.register_container<ContainerTH2Many<float>>("dr_v_eta_second_hits_in_F3",
  360. fv::apply(calc_dr_v_eta, second_hits_in_F3),
  361. "Second Hit in FPIX-L3", params_dr);
  362. //Plots for dz of collections defined above
  363. auto calc_dz_v_eta = lookup_function<pair_vec(vector<HitPair>)>("calc_dz_v_eta");
  364. TH2Params params_dz = {"$\\eta$", 100, -4, 4,
  365. "$\\Delta z$(rad)", 50, -.01, .01};
  366. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_first_hits_in_B1",
  367. fv::apply(calc_dz_v_eta, first_hits_in_B1),
  368. "First Hit in BPIX-L1", params_dz);
  369. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_first_hits_in_B2",
  370. fv::apply(calc_dz_v_eta, first_hits_in_B2),
  371. "First Hit in BPIX-L2", params_dz);
  372. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_first_hits_in_B1_even_ladder",
  373. fv::apply(calc_dz_v_eta, first_hits_in_B1_even_ladder),
  374. "First Hit in BPIX-L1 - Even Ladders", params_dz);
  375. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_first_hits_in_B1_odd_ladder",
  376. fv::apply(calc_dz_v_eta, first_hits_in_B1_odd_ladder),
  377. "First Hit in BPIX-L1 - Odd Ladders", params_dz);
  378. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_first_hits_in_B2_even_ladder",
  379. fv::apply(calc_dz_v_eta, first_hits_in_B2_even_ladder),
  380. "First Hit in BPIX-L2 - Even Ladders", params_dz);
  381. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_first_hits_in_B2_odd_ladder",
  382. fv::apply(calc_dz_v_eta, first_hits_in_B2_odd_ladder),
  383. "First Hit in BPIX-L2 - Odd Ladders", params_dz);
  384. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_second_hits_in_B2",
  385. fv::apply(calc_dz_v_eta, second_hits_in_B2),
  386. "Second Hit in BPIX-L2", params_dz);
  387. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_second_hits_in_B3",
  388. fv::apply(calc_dz_v_eta, second_hits_in_B3),
  389. "Second Hit in BPIX-L3", params_dz);
  390. tds.register_container<ContainerTH2Many<float>>("dz_v_eta_second_hits_in_B4",
  391. fv::apply(calc_dz_v_eta, second_hits_in_B4),
  392. "Second Hit in BPIX-L4", params_dz);
  393. //Plots for drho of collections defined above
  394. auto calc_drho_v_eta = lookup_function<pair_vec(vector<HitPair>)>("calc_drho_v_eta");
  395. TH2Params params_drho = {"$\\eta$", 100, -4, 4,
  396. "$\\Delta \\rho$(cm)", 50, -.01, .01};
  397. tds.register_container<ContainerTH2Many<float>>("drho_v_eta_first_hits_in_F1",
  398. fv::apply(calc_drho_v_eta, first_hits_in_F1),
  399. "First Hit in FPIX-L1", params_drho);
  400. tds.register_container<ContainerTH2Many<float>>("drho_v_eta_first_hits_in_F2",
  401. fv::apply(calc_drho_v_eta, first_hits_in_F2),
  402. "First Hit in FPIX-L2", params_drho);
  403. tds.register_container<ContainerTH2Many<float>>("drho_v_eta_second_hits_in_F2",
  404. fv::apply(calc_drho_v_eta, second_hits_in_F2),
  405. "Second Hit in FPIX-L2", params_drho);
  406. tds.register_container<ContainerTH2Many<float>>("drho_v_eta_second_hits_in_F3",
  407. fv::apply(calc_drho_v_eta, second_hits_in_F3),
  408. "Second Hit in FPIX-L3", params_drho);
  409. }
  410. void setup_sc_residuals(TrackingDataSet& tds){
  411. auto sc_hit_residuals =
  412. func<pair_vec(vector<SuperCluster>,
  413. int, int, int, bool)>("sc_hit_residuals",
  414. FUNC(([](const vector<SuperCluster>& super_clusters,
  415. const int&& det,
  416. const int&& pix_layer,
  417. const int&& hit_number,
  418. const bool&& sel_dRz){
  419. vector<float> residuals;
  420. vector<float> etas;
  421. for(const SuperCluster& super_cluster : super_clusters){ // loop over all supser-clusters
  422. float dRz, dPhi;
  423. unsigned int nSeeds = super_cluster.charge.size();
  424. for(unsigned int seedIdx=0; seedIdx<nSeeds; seedIdx++){ // loop over all seeds associated w/ sc
  425. if(hit_number == 0){
  426. if(super_cluster.lay1[seedIdx] != pix_layer || super_cluster.subDet1[seedIdx] != det) continue;
  427. dRz = super_cluster.dRz1[seedIdx];
  428. dPhi = super_cluster.dPhi1[seedIdx];
  429. }else{
  430. if(super_cluster.lay2[seedIdx] != pix_layer || super_cluster.subDet2[seedIdx] != det) continue;
  431. dRz = super_cluster.dRz2[seedIdx];
  432. dPhi = super_cluster.dPhi2[seedIdx];
  433. }
  434. DEBUG("New Seed, Idx: " << seedIdx << endl
  435. << "\tdRz: " << dRz << endl
  436. << "\tdPhi: " << dPhi << endl);
  437. if(sel_dRz)
  438. residuals.push_back(dRz);
  439. else
  440. residuals.push_back(dPhi);
  441. etas.push_back(pseudorapidityP(super_cluster));
  442. }
  443. }
  444. return std::make_pair(etas,residuals);
  445. })));
  446. auto pix_barrel = constant("PIXEL_BARREL", PIXEL_BARREL);
  447. /* auto pix_endcap = constant<unsigned int>("PIXEL_ENDCAP", PIXEL_ENDCAP); */
  448. auto first_hit = constant("1st", 0);
  449. auto second_hit = constant("2nd", 1);
  450. auto L1 = constant("L1", 1);
  451. auto L2 = constant("L2", 2);
  452. auto L3 = constant("L3", 3);
  453. auto L4 = constant("L4", 4);
  454. auto dRz = constant("dRz", true);
  455. auto dPhi = constant("dPhi", false);
  456. // First hits on inner two bpix layers
  457. auto sc_first_hits_in_B1_dz = fv::tup_apply(sc_hit_residuals,
  458. fv::tuple(super_clusters, pix_barrel, L1, first_hit, dRz), "sc_first_hits_in_B1_dz");
  459. auto sc_first_hits_in_B2_dz = fv::tup_apply(sc_hit_residuals,
  460. fv::tuple(super_clusters, pix_barrel, L2, first_hit, dRz), "sc_first_hits_in_B2_dz");
  461. auto sc_first_hits_in_B3_dz = fv::tup_apply(sc_hit_residuals,
  462. fv::tuple(super_clusters, pix_barrel, L3, first_hit, dRz), "sc_first_hits_in_B3_dz");
  463. auto sc_second_hits_in_B2_dz = fv::tup_apply(sc_hit_residuals,
  464. fv::tuple(super_clusters, pix_barrel, L2, second_hit, dRz), "sc_second_hits_in_B2_dz");
  465. auto sc_second_hits_in_B3_dz = fv::tup_apply(sc_hit_residuals,
  466. fv::tuple(super_clusters, pix_barrel, L3, second_hit, dRz), "sc_second_hits_in_B3_dz");
  467. auto sc_second_hits_in_B4_dz = fv::tup_apply(sc_hit_residuals,
  468. fv::tuple(super_clusters, pix_barrel, L4, second_hit, dRz), "sc_second_hits_in_B4_dz");
  469. auto sc_first_hits_in_B1_dphi = fv::tup_apply(sc_hit_residuals,
  470. fv::tuple(super_clusters, pix_barrel, L1, first_hit, dPhi), "sc_first_hits_in_B1_dphi");
  471. auto sc_first_hits_in_B2_dphi = fv::tup_apply(sc_hit_residuals,
  472. fv::tuple(super_clusters, pix_barrel, L2, first_hit, dPhi), "sc_first_hits_in_B2_dphi");
  473. auto sc_first_hits_in_B3_dphi = fv::tup_apply(sc_hit_residuals,
  474. fv::tuple(super_clusters, pix_barrel, L3, first_hit, dPhi), "sc_first_hits_in_B3_dphi");
  475. auto sc_second_hits_in_B2_dphi = fv::tup_apply(sc_hit_residuals,
  476. fv::tuple(super_clusters, pix_barrel, L2, second_hit, dPhi), "sc_second_hits_in_B2_dphi");
  477. auto sc_second_hits_in_B3_dphi = fv::tup_apply(sc_hit_residuals,
  478. fv::tuple(super_clusters, pix_barrel, L3, second_hit, dPhi), "sc_second_hits_in_B3_dphi");
  479. auto sc_second_hits_in_B4_dphi = fv::tup_apply(sc_hit_residuals,
  480. fv::tuple(super_clusters, pix_barrel, L4, second_hit, dPhi), "sc_second_hits_in_B4_dphi");
  481. TH2Params params_dz = {"$\\eta$", 100, -4, 4,
  482. "$\\delta z$(cm)", 100, -0.7, 0.7};
  483. tds.register_container<ContainerTH2Many<float>>("sc_first_hits_in_B1_dz",
  484. sc_first_hits_in_B1_dz,
  485. "First Hit in BPIX-L1", params_dz);
  486. tds.register_container<ContainerTH2Many<float>>("sc_first_hits_in_B2_dz",
  487. sc_first_hits_in_B2_dz,
  488. "First Hit in BPIX-L2", params_dz);
  489. tds.register_container<ContainerTH2Many<float>>("sc_first_hits_in_B3_dz",
  490. sc_first_hits_in_B3_dz,
  491. "First Hit in BPIX-L3", params_dz);
  492. tds.register_container<ContainerTH2Many<float>>("sc_second_hits_in_B2_dz",
  493. sc_second_hits_in_B2_dz,
  494. "Second Hit in BPIX-L2", params_dz);
  495. tds.register_container<ContainerTH2Many<float>>("sc_second_hits_in_B3_dz",
  496. sc_second_hits_in_B3_dz,
  497. "Second Hit in BPIX-L3", params_dz);
  498. tds.register_container<ContainerTH2Many<float>>("sc_second_hits_in_B4_dz",
  499. sc_second_hits_in_B4_dz,
  500. "Second Hit in BPIX-L4", params_dz);
  501. TH2Params params_dphi = {"$\\eta$", 100, -4, 4,
  502. "$\\delta \\phi$(rad)", 500, -0.7, 0.7};
  503. tds.register_container<ContainerTH2Many<float>>("sc_first_hits_in_B1_dphi",
  504. sc_first_hits_in_B1_dphi,
  505. "First Hit in BPIX-L1", params_dphi);
  506. tds.register_container<ContainerTH2Many<float>>("sc_first_hits_in_B2_dphi",
  507. sc_first_hits_in_B2_dphi,
  508. "First Hit in BPIX-L2", params_dphi);
  509. tds.register_container<ContainerTH2Many<float>>("sc_first_hits_in_B3_dphi",
  510. sc_first_hits_in_B3_dphi,
  511. "First Hit in BPIX-L3", params_dphi);
  512. tds.register_container<ContainerTH2Many<float>>("sc_second_hits_in_B2_dphi",
  513. sc_second_hits_in_B2_dphi,
  514. "Second Hit in BPIX-L2", params_dphi);
  515. tds.register_container<ContainerTH2Many<float>>("sc_second_hits_in_B3_dphi",
  516. sc_second_hits_in_B3_dphi,
  517. "Second Hit in BPIX-L3", params_dphi);
  518. tds.register_container<ContainerTH2Many<float>>("sc_second_hits_in_B4_dphi",
  519. sc_second_hits_in_B4_dphi,
  520. "Second Hit in BPIX-L4", params_dphi);
  521. }
  522. void run_analysis(const vector<fv::util::DataFileDescriptor>& dfds, const string output_filename, bool silent){
  523. gSystem->Load("libfilval.so");
  524. auto replace_suffix = [](const std::string& input, const std::string& new_suffix){
  525. return input.substr(0, input.find_last_of(".")) + new_suffix;
  526. };
  527. string log_filename = replace_suffix(output_filename, ".log");
  528. fv::util::Log::init_logger(log_filename, fv::util::LogPriority::kLogInfo);
  529. TrackingDataSet tds(output_filename, dfds, "trackingNtuple/tree");
  530. /* tds.set_max_events(10); */
  531. register_objects(tds);
  532. setup_hit_pair_functions();
  533. setup_first_hit_pairs(tds);
  534. setup_skipped_layer_hit_pairs(tds);
  535. setup_sc_residuals(tds);
  536. tds.process(silent);
  537. tds.save_all();
  538. }
  539. int main(int argc, char * argv[]){
  540. fv::util::ArgParser args(argc, argv);
  541. bool silent = args.cmdOptionExists("-s");
  542. string output_filename = args.cmdOptionExists("-o") ? args.getCmdOption("-o") : "output.root";
  543. if(args.cmdOptionExists("-F")){
  544. auto file_list = fv::util::read_input_list(args.getCmdOption("-F"));
  545. run_analysis(file_list, output_filename, silent);
  546. }else if(args.cmdOptionExists("-f")){
  547. string input_filename = args.getCmdOption("-f");
  548. string data_label = args.cmdOptionExists("-l") ? args.getCmdOption("-l") : "";
  549. string data_category = args.cmdOptionExists("-c") ? args.getCmdOption("-c") : "";
  550. fv::util::DataFileDescriptor dfd(input_filename, data_label, data_category);
  551. run_analysis(std::vector<fv::util::DataFileDescriptor>({dfd}), output_filename, silent);
  552. }else {
  553. cout << "Usage: ./tracking_validation (-s) {-o output_filename} -F datafiles.txt" << endl;
  554. cout << " ./tracking_validation (-s) {-l DATA_LABEL} {-c DATA_CATEGORY} {-o output_filename} -f treefile.root" << endl;
  555. }
  556. return 0;
  557. }