tracking_validation.cpp 17 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316
  1. #include <iostream>
  2. #include <vector>
  3. #include <utility>
  4. #include <numeric>
  5. #include <limits>
  6. #include <cmath>
  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 "TrackingNtuple.h"
  13. #include "analysis/obj_types.cpp"
  14. using namespace std;
  15. using namespace fv;
  16. using namespace fv::root;
  17. using namespace std;
  18. #define HIT_TYPE_PIXEL 0
  19. #define HIT_TYPE_GLUED 1
  20. #define HIT_TYPE_STRIP 2
  21. #define PIXEL_BARREL 1
  22. #define PIXEL_ENDCAP 2
  23. template<typename A, typename B>
  24. float displacement(const A& a, const B& b){
  25. return std::sqrt(pow(a.x-b.x, 2) +
  26. pow(a.x-b.x, 2) +
  27. pow(a.x-b.x, 2));
  28. }
  29. typedef std::tuple<Seed, PixRecHit, Track, SimHit, SimTrack> MatchedTrack;
  30. typedef std::pair<std::vector<float>,std::vector<float>> pair_vec;
  31. vector<string> seedTypes =
  32. {"initialStepSeeds",
  33. "highPtTripletStepSeeds",
  34. "mixedTripletStepSeeds",
  35. "pixelLessStepSeeds",
  36. "tripletElectronSeeds",
  37. "pixelPairElectronSeeds",
  38. "stripPairElectronSeeds"};
  39. int nSeedTypes = seedTypes.size();
  40. Value<vector<Seed>>* seeds;
  41. Value<vector<PixRecHit>>* pixrec_hits;
  42. Value<vector<Track>>* tracks;
  43. Value<vector<SimHit>>* sim_hits;
  44. Value<vector<SimTrack>>* sim_tracks;
  45. void register_objects(TrackingDataSet& tds){
  46. seeds = register_seeds(tds);
  47. pixrec_hits = register_pixrec_hits(tds);
  48. tracks = register_tracks(tds);
  49. sim_hits = register_sim_hits(tds);
  50. sim_tracks = register_sim_tracks(tds);
  51. }
  52. void setup_matched_tracks(TrackingDataSet& tds){
  53. // Finds pairs of (track,sim_track) where the cycle
  54. // track->seed->rec_hit->sim_hit->sim_track->track
  55. // is closed
  56. auto& find_matched_tracks =
  57. func<vector<MatchedTrack>(vector<Seed>,
  58. vector<PixRecHit>,
  59. vector<Track>,
  60. vector<SimHit>,
  61. vector<SimTrack>)>("find_matched_tracks",
  62. FUNC(([](const vector<Seed>& seeds,
  63. const vector<PixRecHit> pixrec_hits,
  64. const vector<Track>& tracks,
  65. const vector<SimHit> sim_hits,
  66. const vector<SimTrack> sim_tracks){
  67. INFO("New event");
  68. INFO(boost::format("Number of tracks is %d, number of seeds is %d, number of hits %d, number of simhits %d")
  69. % tracks.size() % seeds.size() % pixrec_hits.size() % sim_hits.size());
  70. vector<MatchedTrack> matched_tracks;
  71. for(const Track &track : tracks){
  72. const Seed& seed = seeds[track.seedIdx];
  73. if( seed.algoOriginal < 0 || seed.algoOriginal >= nSeedTypes) continue;
  74. INFO(boost::format(" seed indx= %d algorithm= %d type %s\n")
  75. % track.seedIdx % seed.algoOriginal % seedTypes[seed.algoOriginal]);
  76. for(auto tup : boost::combine(seed.hitIdx, seed.hitType)){ //for hits in track's seed
  77. int hitIdx, hitType;
  78. boost::tie(hitIdx, hitType) = tup;
  79. if(hitType != HIT_TYPE_PIXEL) continue; // take only pixel hits for now
  80. const PixRecHit &rec_hit = pixrec_hits[hitIdx];
  81. TVector3 recoHitPoint;
  82. recoHitPoint.SetXYZ( rec_hit.x, rec_hit.y, rec_hit.z);
  83. float Rreco = recoHitPoint.Perp();
  84. INFO(boost::format(" hit= %d is at radius= %f\n") % hitIdx % Rreco);
  85. for(const int& simHitIdx : rec_hit.simHitIdx){ // for sim-hits matching rec-hit
  86. const SimHit &sim_hit = sim_hits[simHitIdx];
  87. const SimTrack &sim_track = sim_tracks[sim_hit.simTrkIdx];
  88. for(const int &trkIdx : sim_track.trkIdx){ // for tracks matched to this sim_track
  89. if(trkIdx == track.idx){ // sim_track matches with our original track
  90. matched_tracks.push_back({seed, rec_hit, track, sim_hit, sim_track});
  91. break;
  92. }
  93. }
  94. }
  95. }
  96. }
  97. return matched_tracks;
  98. })));
  99. fv::apply(find_matched_tracks, fv::tuple(seeds, pixrec_hits, tracks, sim_hits, sim_tracks), "matched_tracks");
  100. }
  101. bool in_det(const MatchedTrack &mt, int &&det, int &&layer){
  102. auto& hit = std::get<PixRecHit>(mt);
  103. return hit.det == det && hit.lay == layer;
  104. };
  105. void setup_residuals(TrackingDataSet &tds){
  106. auto matched_tracks = lookup<vector<MatchedTrack>>("matched_tracks");
  107. // Break matched_tracks into catagories based on Barrel/Endcap and layer
  108. auto matched_tracks_B1 = filter(func<bool(MatchedTrack)>("matched_tracks_B1",
  109. FUNC(([](const MatchedTrack& matched_track){
  110. return in_det(matched_track, PIXEL_BARREL, 1);
  111. }))), matched_tracks);
  112. auto matched_tracks_B2 = filter(func<bool(MatchedTrack)>("matched_tracks_B2",
  113. FUNC(([](const MatchedTrack& matched_track){
  114. return in_det(matched_track, PIXEL_BARREL, 2);
  115. }))), matched_tracks);
  116. auto matched_tracks_B3 = filter(func<bool(MatchedTrack)>("matched_tracks_B3",
  117. FUNC(([](const MatchedTrack& matched_track){
  118. return in_det(matched_track, PIXEL_BARREL, 3);
  119. }))), matched_tracks);
  120. auto matched_tracks_B4 = filter(func<bool(MatchedTrack)>("matched_tracks_B4",
  121. FUNC(([](const MatchedTrack& matched_track){
  122. return in_det(matched_track, PIXEL_BARREL, 4);
  123. }))), matched_tracks);
  124. auto matched_tracks_F1 = filter(func<bool(MatchedTrack)>("matched_tracks_F1",
  125. FUNC(([](const MatchedTrack& matched_track){
  126. return in_det(matched_track, PIXEL_ENDCAP, 1);
  127. }))), matched_tracks);
  128. auto matched_tracks_F2 = filter(func<bool(MatchedTrack)>("matched_tracks_F2",
  129. FUNC(([](const MatchedTrack& matched_track){
  130. return in_det(matched_track, PIXEL_ENDCAP, 2);
  131. }))), matched_tracks);
  132. auto matched_tracks_F3 = filter(func<bool(MatchedTrack)>("matched_tracks_F3",
  133. FUNC(([](const MatchedTrack& matched_track){
  134. return in_det(matched_track, PIXEL_ENDCAP, 3);
  135. }))), matched_tracks);
  136. auto& calc_residuals_v_eta = func<pair_vec(vector<MatchedTrack>)>("matched_track_hit_residuals",
  137. FUNC(([](const vector<MatchedTrack>& matched_tracks){
  138. vector<float> residuals;
  139. vector<float> etas;
  140. for(auto matched_track : matched_tracks){
  141. auto& pixrec_hit = std::get<PixRecHit>(matched_track);
  142. auto& sim_hit = std::get<SimHit>(matched_track);
  143. auto& sim_track = std::get<SimTrack>(matched_track);
  144. residuals.push_back(displacement(pixrec_hit, sim_hit));
  145. etas.push_back(sim_track.eta);
  146. }
  147. return std::make_pair(etas, residuals);
  148. })));
  149. tds.register_container<ContainerTH2Many<float>>("matched_track_residuals_v_eta_all", "Matched Track Residuals - All",
  150. fv::apply(calc_residuals_v_eta, fv::tuple(matched_tracks)),
  151. 100, -4, 4, 100, 0, .01,"Eta","Residuals(cm)");
  152. tds.register_container<ContainerTH2Many<float>>("matched_track_residuals_v_eta_B1", "Matched Track Residuals - B1",
  153. fv::apply(calc_residuals_v_eta, fv::tuple(matched_tracks_B1)),
  154. 100, -4, 4, 100, 0, .01,"Eta","Residuals(cm)");
  155. tds.register_container<ContainerTH2Many<float>>("matched_track_residuals_v_eta_B2", "Matched Track Residuals - B2",
  156. fv::apply(calc_residuals_v_eta, fv::tuple(matched_tracks_B2)),
  157. 100, -4, 4, 100, 0, .01,"Eta","Residuals(cm)");
  158. tds.register_container<ContainerTH2Many<float>>("matched_track_residuals_v_eta_B3", "Matched Track Residuals - B3",
  159. fv::apply(calc_residuals_v_eta, fv::tuple(matched_tracks_B3)),
  160. 100, -4, 4, 100, 0, .01,"Eta","Residuals(cm)");
  161. tds.register_container<ContainerTH2Many<float>>("matched_track_residuals_v_eta_B4", "Matched Track Residuals - B4",
  162. fv::apply(calc_residuals_v_eta, fv::tuple(matched_tracks_B4)),
  163. 100, -4, 4, 100, 0, .01,"Eta","Residuals(cm)");
  164. tds.register_container<ContainerTH2Many<float>>("matched_track_residuals_v_eta_F1", "Matched Track Residuals - F1",
  165. fv::apply(calc_residuals_v_eta, fv::tuple(matched_tracks_F1)),
  166. 100, -4, 4, 100, 0, .01,"Eta","Residuals(cm)");
  167. tds.register_container<ContainerTH2Many<float>>("matched_track_residuals_v_eta_F2", "Matched Track Residuals - F2",
  168. fv::apply(calc_residuals_v_eta, fv::tuple(matched_tracks_F2)),
  169. 100, -4, 4, 100, 0, .01,"Eta","Residuals(cm)");
  170. tds.register_container<ContainerTH2Many<float>>("matched_track_residuals_v_eta_F3", "Matched Track Residuals - F3",
  171. fv::apply(calc_residuals_v_eta, fv::tuple(matched_tracks_F3)),
  172. 100, -4, 4, 100, 0, .01,"Eta","Residuals(cm)");
  173. auto& calc_dphi_v_eta = func<pair_vec(vector<MatchedTrack>)>("matched_track_hit_dphis",
  174. FUNC(([](const vector<MatchedTrack>& matched_tracks){
  175. vector<float> dphis;
  176. vector<float> etas;
  177. for(auto matched_track : matched_tracks){
  178. auto& pixrec_hit = std::get<PixRecHit>(matched_track);
  179. auto& sim_hit = std::get<SimHit>(matched_track);
  180. auto& sim_track = std::get<SimTrack>(matched_track);
  181. dphis.push_back(atan2(sim_hit.x, sim_hit.y) - atan2(pixrec_hit.x, pixrec_hit.y));
  182. etas.push_back(sim_track.eta);
  183. }
  184. return std::make_pair(etas, dphis);
  185. })));
  186. tds.register_container<ContainerTH2Many<float>>("matched_track_dphis_v_eta_all", "Matched Track dphis - All",
  187. fv::apply(calc_dphi_v_eta, fv::tuple(matched_tracks)),
  188. 100, -4, 4, 100, -.002, .002,"Eta","Residuals(rad)");
  189. tds.register_container<ContainerTH2Many<float>>("matched_track_dphis_v_eta_B1", "Matched Track dphis - B1",
  190. fv::apply(calc_dphi_v_eta, fv::tuple(matched_tracks_B1)),
  191. 100, -4, 4, 100, -.002, .002,"Eta","dphis(rad)");
  192. tds.register_container<ContainerTH2Many<float>>("matched_track_dphis_v_eta_B2", "Matched Track dphis - B2",
  193. fv::apply(calc_dphi_v_eta, fv::tuple(matched_tracks_B2)),
  194. 100, -4, 4, 100, -.002, .002,"Eta","dphis(rad)");
  195. tds.register_container<ContainerTH2Many<float>>("matched_track_dphis_v_eta_B3", "Matched Track dphis - B3",
  196. fv::apply(calc_dphi_v_eta, fv::tuple(matched_tracks_B3)),
  197. 100, -4, 4, 100, -.002, .002,"Eta","dphis(rad)");
  198. tds.register_container<ContainerTH2Many<float>>("matched_track_dphis_v_eta_B4", "Matched Track dphis - B4",
  199. fv::apply(calc_dphi_v_eta, fv::tuple(matched_tracks_B4)),
  200. 100, -4, 4, 100, -.002, .002,"Eta","dphis(rad)");
  201. tds.register_container<ContainerTH2Many<float>>("matched_track_dphis_v_eta_F1", "Matched Track dphis - F1",
  202. fv::apply(calc_dphi_v_eta, fv::tuple(matched_tracks_F1)),
  203. 100, -4, 4, 100, -.002, .002,"Eta","dphis(rad)");
  204. tds.register_container<ContainerTH2Many<float>>("matched_track_dphis_v_eta_F2", "Matched Track dphis - F2",
  205. fv::apply(calc_dphi_v_eta, fv::tuple(matched_tracks_F2)),
  206. 100, -4, 4, 100, -.002, .002,"Eta","dphis(rad)");
  207. tds.register_container<ContainerTH2Many<float>>("matched_track_dphis_v_eta_F3", "Matched Track dphis - F3",
  208. fv::apply(calc_dphi_v_eta, fv::tuple(matched_tracks_F3)),
  209. 100, -4, 4, 100, -.002, .002,"Eta","dphis(rad)");
  210. auto& calc_dz_v_eta = func<pair_vec(vector<MatchedTrack>)>("matched_track_hit_dphis",
  211. FUNC(([](const vector<MatchedTrack>& matched_tracks){
  212. vector<float> dzs;
  213. vector<float> etas;
  214. for(auto matched_track : matched_tracks){
  215. auto& pixrec_hit = std::get<PixRecHit>(matched_track);
  216. auto& sim_hit = std::get<SimHit>(matched_track);
  217. auto& sim_track = std::get<SimTrack>(matched_track);
  218. dzs.push_back(sim_hit.z - pixrec_hit.z);
  219. etas.push_back(sim_track.eta);
  220. }
  221. return std::make_pair(etas, dzs);
  222. })));
  223. tds.register_container<ContainerTH2Many<float>>("matched_track_dzs_v_eta_all", "Matched Track dz - All",
  224. fv::apply(calc_dz_v_eta, fv::tuple(matched_tracks)),
  225. 100, -4, 4, 100, -.002, .002,"Eta","dz(cm)");
  226. tds.register_container<ContainerTH2Many<float>>("matched_track_dzs_v_eta_B1", "Matched Track dz - B1",
  227. fv::apply(calc_dz_v_eta, fv::tuple(matched_tracks_B1)),
  228. 100, -4, 4, 100, -.002, .002,"Eta","dz(cm)");
  229. tds.register_container<ContainerTH2Many<float>>("matched_track_dzs_v_eta_B2", "Matched Track dz - B2",
  230. fv::apply(calc_dz_v_eta, fv::tuple(matched_tracks_B2)),
  231. 100, -4, 4, 100, -.002, .002,"Eta","dz(cm)");
  232. tds.register_container<ContainerTH2Many<float>>("matched_track_dzs_v_eta_B3", "Matched Track dz - B3",
  233. fv::apply(calc_dz_v_eta, fv::tuple(matched_tracks_B3)),
  234. 100, -4, 4, 100, -.002, .002,"Eta","dz(cm)");
  235. tds.register_container<ContainerTH2Many<float>>("matched_track_dzs_v_eta_B4", "Matched Track dz - B4",
  236. fv::apply(calc_dz_v_eta, fv::tuple(matched_tracks_B4)),
  237. 100, -4, 4, 100, -.002, .002,"Eta","dz(cm)");
  238. tds.register_container<ContainerTH2Many<float>>("matched_track_dzs_v_eta_F1", "Matched Track dz - F1",
  239. fv::apply(calc_dz_v_eta, fv::tuple(matched_tracks_F1)),
  240. 100, -4, 4, 100, -.002, .002,"Eta","dz(cm)");
  241. tds.register_container<ContainerTH2Many<float>>("matched_track_dzs_v_eta_F2", "Matched Track dz - F2",
  242. fv::apply(calc_dz_v_eta, fv::tuple(matched_tracks_F2)),
  243. 100, -4, 4, 100, -.002, .002,"Eta","dz(cm)");
  244. tds.register_container<ContainerTH2Many<float>>("matched_track_dzs_v_eta_F3", "Matched Track dz - F3",
  245. fv::apply(calc_dz_v_eta, fv::tuple(matched_tracks_F3)),
  246. 100, -4, 4, 100, -.002, .002,"Eta","dz(cm)");
  247. }
  248. void run_analysis(const vector<fv::util::DataFileDescriptor>& dfds, const string output_filename, bool silent){
  249. gSystem->Load("libfilval.so");
  250. auto replace_suffix = [](const std::string& input, const std::string& new_suffix){
  251. return input.substr(0, input.find_last_of(".")) + new_suffix;
  252. };
  253. string log_filename = replace_suffix(output_filename, ".log");
  254. fv::util::Log::init_logger(log_filename, fv::util::LogPriority::kLogDebug);
  255. TrackingDataSet tds(output_filename, dfds, "trackingNtuple/tree");
  256. /* tds.set_max_events(10); */
  257. register_objects(tds);
  258. setup_matched_tracks(tds);
  259. setup_residuals(tds);
  260. tds.process(silent);
  261. tds.save_all();
  262. }
  263. int main(int argc, char * argv[])
  264. {
  265. fv::util::ArgParser args(argc, argv);
  266. bool silent = args.cmdOptionExists("-s");
  267. string output_filename = args.cmdOptionExists("-o") ? args.getCmdOption("-o") : "output.root";
  268. if(args.cmdOptionExists("-F")){
  269. auto file_list = fv::util::read_input_list(args.getCmdOption("-F"));
  270. run_analysis(file_list, output_filename, silent);
  271. }else if(args.cmdOptionExists("-f")){
  272. string input_filename = args.getCmdOption("-f");
  273. string data_label = args.cmdOptionExists("-l") ? args.getCmdOption("-l") : "";
  274. fv::util::DataFileDescriptor dfd(input_filename, data_label);
  275. run_analysis(std::vector<fv::util::DataFileDescriptor>({dfd}), output_filename, silent);
  276. }else {
  277. cout << "Usage: ./tracking_validation (-s) {-o output_filename} -F datafiles.txt" << endl;
  278. cout << " ./tracking_validation (-s) {-l DATA_LABEL} {-o output_filename} -f treefile.root" << endl;
  279. }
  280. }