tracking_validation.cpp 36 KB

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