CMS 3D CMS Logo

TICLDumper.cc
Go to the documentation of this file.
1 // Original Authors: Philipp Zehetner, Wahid Redjeb
2 
3 #include "TTree.h"
4 #include "TFile.h"
5 
6 #include <iostream>
7 #include <fstream>
8 #include <sstream>
9 #include <variant>
10 
11 #include <memory> // unique_ptr
23 
36 
42 
45 
49 
53 
55 
56 // TFileService
59 
60 class TICLDumper : public edm::one::EDAnalyzer<edm::one::WatchRuns, edm::one::SharedResources> {
61 public:
62  explicit TICLDumper(const edm::ParameterSet&);
63  ~TICLDumper() override;
64  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
66  typedef std::vector<double> Vec;
67 
68 private:
69  void beginJob() override;
70  void beginRun(const edm::Run&, const edm::EventSetup&) override;
71 
72  void initialize(const HGCalDDDConstants* hgcons,
73  const hgcal::RecHitTools rhtools,
74  const edm::ESHandle<MagneticField> bfieldH,
75  const edm::ESHandle<Propagator> propH);
76  void buildLayers();
77 
78  void analyze(const edm::Event&, const edm::EventSetup&) override;
79  void endRun(edm::Run const& iEvent, edm::EventSetup const&) override{};
80  void endJob() override;
81 
82  // Define Tokens
125 
134  std::unique_ptr<GeomDet> firstDisk_[2];
135  std::unique_ptr<GeomDet> interfaceDisk_[2];
138  bool saveLCs_;
147 
148  // Output tree
149  TTree* tree_;
150 
151  void clearVariables();
152 
153  // Variables for branches
154  unsigned int ev_event_;
155  unsigned int ntracksters_;
156  unsigned int nclusters_;
157  unsigned int stsSC_ntracksters_;
158  unsigned int stsCP_ntracksters_;
161 
162  std::vector<float> trackster_time;
163  std::vector<float> trackster_timeError;
164  std::vector<float> trackster_regressed_energy;
165  std::vector<float> trackster_raw_energy;
166  std::vector<float> trackster_raw_em_energy;
167  std::vector<float> trackster_raw_pt;
168  std::vector<float> trackster_raw_em_pt;
169  std::vector<float> trackster_barycenter_x;
170  std::vector<float> trackster_barycenter_y;
171  std::vector<float> trackster_barycenter_z;
172  std::vector<float> trackster_EV1;
173  std::vector<float> trackster_EV2;
174  std::vector<float> trackster_EV3;
175  std::vector<float> trackster_eVector0_x;
176  std::vector<float> trackster_eVector0_y;
177  std::vector<float> trackster_eVector0_z;
178  std::vector<float> trackster_sigmaPCA1;
179  std::vector<float> trackster_sigmaPCA2;
180  std::vector<float> trackster_sigmaPCA3;
181  std::vector<float> trackster_barycenter_eta;
182  std::vector<float> trackster_barycenter_phi;
183  std::vector<std::vector<float>> trackster_id_probabilities;
184  std::vector<std::vector<uint32_t>> trackster_vertices_indexes;
185  std::vector<std::vector<float>> trackster_vertices_x;
186  std::vector<std::vector<float>> trackster_vertices_y;
187  std::vector<std::vector<float>> trackster_vertices_z;
188  std::vector<std::vector<float>> trackster_vertices_time;
189  std::vector<std::vector<float>> trackster_vertices_timeErr;
190  std::vector<std::vector<float>> trackster_vertices_energy;
191  std::vector<std::vector<float>> trackster_vertices_correctedEnergy;
192  std::vector<std::vector<float>> trackster_vertices_correctedEnergyUncertainty;
193  std::vector<std::vector<float>> trackster_vertices_multiplicity;
194 
195  std::vector<float> stsSC_trackster_time;
196  std::vector<float> stsSC_trackster_timeBoundary;
197  std::vector<float> stsSC_trackster_timeError;
199  std::vector<float> stsSC_trackster_regressed_pt;
200  std::vector<float> stsSC_trackster_raw_energy;
201  std::vector<float> stsSC_trackster_raw_em_energy;
202  std::vector<float> stsSC_trackster_raw_pt;
203  std::vector<float> stsSC_trackster_raw_em_pt;
204  std::vector<float> stsSC_trackster_barycenter_x;
205  std::vector<float> stsSC_trackster_barycenter_y;
206  std::vector<float> stsSC_trackster_barycenter_z;
207  std::vector<float> stsSC_trackster_barycenter_eta;
208  std::vector<float> stsSC_trackster_barycenter_phi;
209  std::vector<float> stsSC_trackster_EV1;
210  std::vector<float> stsSC_trackster_EV2;
211  std::vector<float> stsSC_trackster_EV3;
212  std::vector<float> stsSC_trackster_eVector0_x;
213  std::vector<float> stsSC_trackster_eVector0_y;
214  std::vector<float> stsSC_trackster_eVector0_z;
215  std::vector<float> stsSC_trackster_sigmaPCA1;
216  std::vector<float> stsSC_trackster_sigmaPCA2;
217  std::vector<float> stsSC_trackster_sigmaPCA3;
218  std::vector<int> stsSC_pdgID;
219  std::vector<int> stsSC_trackIdx;
220  std::vector<float> stsSC_trackTime;
221  std::vector<float> stsSC_boundaryX;
222  std::vector<float> stsSC_boundaryY;
223  std::vector<float> stsSC_boundaryZ;
224  std::vector<float> stsSC_boundaryEta;
225  std::vector<float> stsSC_boundaryPhi;
226  std::vector<float> stsSC_boundaryPx;
227  std::vector<float> stsSC_boundaryPy;
228  std::vector<float> stsSC_boundaryPz;
229  std::vector<float> stsSC_track_boundaryX;
230  std::vector<float> stsSC_track_boundaryY;
231  std::vector<float> stsSC_track_boundaryZ;
232  std::vector<float> stsSC_track_boundaryEta;
233  std::vector<float> stsSC_track_boundaryPhi;
234  std::vector<float> stsSC_track_boundaryPx;
235  std::vector<float> stsSC_track_boundaryPy;
236  std::vector<float> stsSC_track_boundaryPz;
237  std::vector<std::vector<float>> stsSC_trackster_id_probabilities;
238  std::vector<std::vector<uint32_t>> stsSC_trackster_vertices_indexes;
239  std::vector<std::vector<float>> stsSC_trackster_vertices_x;
240  std::vector<std::vector<float>> stsSC_trackster_vertices_y;
241  std::vector<std::vector<float>> stsSC_trackster_vertices_z;
242  std::vector<std::vector<float>> stsSC_trackster_vertices_time;
243  std::vector<std::vector<float>> stsSC_trackster_vertices_timeErr;
244  std::vector<std::vector<float>> stsSC_trackster_vertices_energy;
245  std::vector<std::vector<float>> stsSC_trackster_vertices_correctedEnergy;
247  std::vector<std::vector<float>> stsSC_trackster_vertices_multiplicity;
248  std::vector<float> stsCP_trackster_time;
249  std::vector<float> stsCP_trackster_timeBoundary;
250  std::vector<float> stsCP_trackster_timeError;
252  std::vector<float> stsCP_trackster_regressed_pt;
253  std::vector<float> stsCP_trackster_raw_energy;
254  std::vector<float> stsCP_trackster_raw_em_energy;
255  std::vector<float> stsCP_trackster_raw_pt;
256  std::vector<float> stsCP_trackster_raw_em_pt;
257  std::vector<float> stsCP_trackster_barycenter_x;
258  std::vector<float> stsCP_trackster_barycenter_y;
259  std::vector<float> stsCP_trackster_barycenter_z;
260  std::vector<float> stsCP_trackster_barycenter_eta;
261  std::vector<float> stsCP_trackster_barycenter_phi;
262  std::vector<float> stsCP_trackster_EV1;
263  std::vector<float> stsCP_trackster_EV2;
264  std::vector<float> stsCP_trackster_EV3;
265  std::vector<float> stsCP_trackster_eVector0_x;
266  std::vector<float> stsCP_trackster_eVector0_y;
267  std::vector<float> stsCP_trackster_eVector0_z;
268  std::vector<float> stsCP_trackster_sigmaPCA1;
269  std::vector<float> stsCP_trackster_sigmaPCA2;
270  std::vector<float> stsCP_trackster_sigmaPCA3;
271  std::vector<int> stsCP_pdgID;
272  std::vector<int> stsCP_trackIdx;
273  std::vector<float> stsCP_trackTime;
274  std::vector<float> stsCP_boundaryX;
275  std::vector<float> stsCP_boundaryY;
276  std::vector<float> stsCP_boundaryZ;
277  std::vector<float> stsCP_boundaryEta;
278  std::vector<float> stsCP_boundaryPhi;
279  std::vector<float> stsCP_boundaryPx;
280  std::vector<float> stsCP_boundaryPy;
281  std::vector<float> stsCP_boundaryPz;
282  std::vector<float> stsCP_track_boundaryX;
283  std::vector<float> stsCP_track_boundaryY;
284  std::vector<float> stsCP_track_boundaryZ;
285  std::vector<float> stsCP_track_boundaryEta;
286  std::vector<float> stsCP_track_boundaryPhi;
287  std::vector<float> stsCP_track_boundaryPx;
288  std::vector<float> stsCP_track_boundaryPy;
289  std::vector<float> stsCP_track_boundaryPz;
290  std::vector<std::vector<float>> stsCP_trackster_id_probabilities;
291  std::vector<std::vector<uint32_t>> stsCP_trackster_vertices_indexes;
292  std::vector<std::vector<float>> stsCP_trackster_vertices_x;
293  std::vector<std::vector<float>> stsCP_trackster_vertices_y;
294  std::vector<std::vector<float>> stsCP_trackster_vertices_z;
295  std::vector<std::vector<float>> stsCP_trackster_vertices_time;
296  std::vector<std::vector<float>> stsCP_trackster_vertices_timeErr;
297  std::vector<std::vector<float>> stsCP_trackster_vertices_energy;
298  std::vector<std::vector<float>> stsCP_trackster_vertices_correctedEnergy;
300  std::vector<std::vector<float>> stsCP_trackster_vertices_multiplicity;
301 
302  std::vector<float> simTICLCandidate_raw_energy;
304  std::vector<std::vector<int>> simTICLCandidate_simTracksterCPIndex;
305  std::vector<float> simTICLCandidate_boundaryX;
306  std::vector<float> simTICLCandidate_boundaryY;
307  std::vector<float> simTICLCandidate_boundaryZ;
308  std::vector<float> simTICLCandidate_boundaryPx;
309  std::vector<float> simTICLCandidate_boundaryPy;
310  std::vector<float> simTICLCandidate_boundaryPz;
312  std::vector<float> simTICLCandidate_time;
313  std::vector<int> simTICLCandidate_pdgId;
314  std::vector<int> simTICLCandidate_charge;
316 
317  // from TICLCandidate, product of linking
318  size_t nCandidates;
319  std::vector<int> candidate_charge;
320  std::vector<int> candidate_pdgId;
321  std::vector<float> candidate_energy;
322  std::vector<float> candidate_raw_energy;
323  std::vector<double> candidate_px;
324  std::vector<double> candidate_py;
325  std::vector<double> candidate_pz;
326  std::vector<float> candidate_time;
327  std::vector<float> candidate_time_err;
328  std::vector<std::vector<float>> candidate_id_probabilities;
329  std::vector<std::vector<uint32_t>> tracksters_in_candidate;
330  std::vector<int> track_in_candidate;
331 
332  // merged tracksters
334  std::vector<float> tracksters_merged_time;
335  std::vector<float> tracksters_merged_timeError;
337  std::vector<float> tracksters_merged_raw_energy;
339  std::vector<float> tracksters_merged_raw_pt;
340  std::vector<float> tracksters_merged_raw_em_pt;
341  std::vector<float> tracksters_merged_barycenter_x;
342  std::vector<float> tracksters_merged_barycenter_y;
343  std::vector<float> tracksters_merged_barycenter_z;
346  std::vector<float> tracksters_merged_EV1;
347  std::vector<float> tracksters_merged_EV2;
348  std::vector<float> tracksters_merged_EV3;
349  std::vector<float> tracksters_merged_eVector0_x;
350  std::vector<float> tracksters_merged_eVector0_y;
351  std::vector<float> tracksters_merged_eVector0_z;
352  std::vector<float> tracksters_merged_sigmaPCA1;
353  std::vector<float> tracksters_merged_sigmaPCA2;
354  std::vector<float> tracksters_merged_sigmaPCA3;
355  std::vector<std::vector<uint32_t>> tracksters_merged_vertices_indexes;
356  std::vector<std::vector<float>> tracksters_merged_vertices_x;
357  std::vector<std::vector<float>> tracksters_merged_vertices_y;
358  std::vector<std::vector<float>> tracksters_merged_vertices_z;
359  std::vector<std::vector<float>> tracksters_merged_vertices_time;
360  std::vector<std::vector<float>> tracksters_merged_vertices_timeErr;
361  std::vector<std::vector<float>> tracksters_merged_vertices_energy;
362  std::vector<std::vector<float>> tracksters_merged_vertices_correctedEnergy;
364  std::vector<std::vector<float>> tracksters_merged_vertices_multiplicity;
365  std::vector<std::vector<float>> tracksters_merged_id_probabilities;
366 
367  // associations
368  std::vector<std::vector<uint32_t>> trackstersCLUE3D_recoToSim_SC;
369  std::vector<std::vector<float>> trackstersCLUE3D_recoToSim_SC_score;
370  std::vector<std::vector<float>> trackstersCLUE3D_recoToSim_SC_sharedE;
371  std::vector<std::vector<uint32_t>> trackstersCLUE3D_simToReco_SC;
372  std::vector<std::vector<float>> trackstersCLUE3D_simToReco_SC_score;
373  std::vector<std::vector<float>> trackstersCLUE3D_simToReco_SC_sharedE;
374 
375  std::vector<std::vector<uint32_t>> trackstersCLUE3D_recoToSim_CP;
376  std::vector<std::vector<float>> trackstersCLUE3D_recoToSim_CP_score;
377  std::vector<std::vector<float>> trackstersCLUE3D_recoToSim_CP_sharedE;
378  std::vector<std::vector<uint32_t>> trackstersCLUE3D_simToReco_CP;
379  std::vector<std::vector<float>> trackstersCLUE3D_simToReco_CP_score;
380  std::vector<std::vector<float>> trackstersCLUE3D_simToReco_CP_sharedE;
381 
382  std::vector<std::vector<uint32_t>> MergeTracksters_recoToSim_SC;
383  std::vector<std::vector<float>> MergeTracksters_recoToSim_SC_score;
384  std::vector<std::vector<float>> MergeTracksters_recoToSim_SC_sharedE;
385  std::vector<std::vector<uint32_t>> MergeTracksters_simToReco_SC;
386  std::vector<std::vector<float>> MergeTracksters_simToReco_SC_score;
387  std::vector<std::vector<float>> MergeTracksters_simToReco_SC_sharedE;
388 
389  std::vector<std::vector<uint32_t>> MergeTracksters_recoToSim_CP;
390  std::vector<std::vector<float>> MergeTracksters_recoToSim_CP_score;
391  std::vector<std::vector<float>> MergeTracksters_recoToSim_CP_sharedE;
392  std::vector<std::vector<uint32_t>> MergeTracksters_simToReco_CP;
393  std::vector<std::vector<float>> MergeTracksters_simToReco_CP_score;
394  std::vector<std::vector<float>> MergeTracksters_simToReco_CP_sharedE;
395 
396  std::vector<std::vector<uint32_t>> MergeTracksters_recoToSim_PU;
397  std::vector<std::vector<float>> MergeTracksters_recoToSim_PU_score;
398  std::vector<std::vector<float>> MergeTracksters_recoToSim_PU_sharedE;
399  std::vector<std::vector<uint32_t>> MergeTracksters_simToReco_PU;
400  std::vector<std::vector<float>> MergeTracksters_simToReco_PU_score;
401  std::vector<std::vector<float>> MergeTracksters_simToReco_PU_sharedE;
402 
403  std::vector<uint32_t> cluster_seedID;
404  std::vector<float> cluster_energy;
405  std::vector<float> cluster_correctedEnergy;
407  std::vector<float> cluster_position_x;
408  std::vector<float> cluster_position_y;
409  std::vector<float> cluster_position_z;
410  std::vector<float> cluster_position_eta;
411  std::vector<float> cluster_position_phi;
412  std::vector<unsigned int> cluster_layer_id;
413  std::vector<int> cluster_type;
414  std::vector<float> cluster_time;
415  std::vector<float> cluster_timeErr;
416  std::vector<uint32_t> cluster_number_of_hits;
417 
418  std::vector<unsigned int> track_id;
419  std::vector<float> track_hgcal_x;
420  std::vector<float> track_hgcal_y;
421  std::vector<float> track_hgcal_z;
422  std::vector<float> track_hgcal_px;
423  std::vector<float> track_hgcal_py;
424  std::vector<float> track_hgcal_pz;
425  std::vector<float> track_hgcal_eta;
426  std::vector<float> track_hgcal_phi;
427  std::vector<float> track_hgcal_pt;
428  std::vector<float> track_pt;
429  std::vector<int> track_quality;
430  std::vector<int> track_missing_outer_hits;
431  std::vector<int> track_missing_inner_hits;
432  std::vector<int> track_charge;
433  std::vector<double> track_time;
434  std::vector<float> track_time_quality;
435  std::vector<float> track_time_err;
436  std::vector<float> track_beta;
437  std::vector<float> track_time_mtd;
438  std::vector<float> track_time_mtd_err;
439  std::vector<GlobalPoint> track_pos_mtd;
440  std::vector<int> track_nhits;
441  std::vector<int> track_isMuon;
442  std::vector<int> track_isTrackerMuon;
443 
451  TTree* tracks_tree_;
453 };
454 
456  // event info
457  ntracksters_ = 0;
458  nclusters_ = 0;
459 
460  trackster_time.clear();
461  trackster_timeError.clear();
463  trackster_raw_energy.clear();
464  trackster_raw_em_energy.clear();
465  trackster_raw_pt.clear();
466  trackster_raw_em_pt.clear();
467  trackster_barycenter_x.clear();
468  trackster_barycenter_y.clear();
469  trackster_barycenter_z.clear();
470  trackster_EV1.clear();
471  trackster_EV2.clear();
472  trackster_EV3.clear();
473  trackster_eVector0_x.clear();
474  trackster_eVector0_y.clear();
475  trackster_eVector0_z.clear();
476  trackster_sigmaPCA1.clear();
477  trackster_sigmaPCA2.clear();
478  trackster_sigmaPCA3.clear();
479  trackster_barycenter_eta.clear();
480  trackster_barycenter_phi.clear();
483  trackster_vertices_x.clear();
484  trackster_vertices_y.clear();
485  trackster_vertices_z.clear();
486  trackster_vertices_time.clear();
492 
493  stsSC_trackster_time.clear();
500  stsSC_trackster_raw_pt.clear();
505  stsSC_trackster_EV1.clear();
506  stsSC_trackster_EV2.clear();
507  stsSC_trackster_EV3.clear();
516  stsSC_pdgID.clear();
517  stsSC_trackIdx.clear();
518  stsSC_trackTime.clear();
519  stsSC_boundaryX.clear();
520  stsSC_boundaryY.clear();
521  stsSC_boundaryZ.clear();
522  stsSC_boundaryEta.clear();
523  stsSC_boundaryPhi.clear();
524  stsSC_boundaryPx.clear();
525  stsSC_boundaryPy.clear();
526  stsSC_boundaryPz.clear();
527  stsSC_track_boundaryX.clear();
528  stsSC_track_boundaryY.clear();
529  stsSC_track_boundaryZ.clear();
530  stsSC_track_boundaryEta.clear();
531  stsSC_track_boundaryPhi.clear();
532  stsSC_track_boundaryPx.clear();
533  stsSC_track_boundaryPy.clear();
534  stsSC_track_boundaryPz.clear();
546 
547  stsCP_trackster_time.clear();
554  stsCP_trackster_raw_pt.clear();
564  stsCP_pdgID.clear();
565  stsCP_trackIdx.clear();
566  stsCP_trackTime.clear();
567  stsCP_boundaryX.clear();
568  stsCP_boundaryY.clear();
569  stsCP_boundaryZ.clear();
570  stsCP_boundaryEta.clear();
571  stsCP_boundaryPhi.clear();
572  stsCP_boundaryPx.clear();
573  stsCP_boundaryPy.clear();
574  stsCP_boundaryPz.clear();
575  stsCP_track_boundaryX.clear();
576  stsCP_track_boundaryY.clear();
577  stsCP_track_boundaryZ.clear();
578  stsCP_track_boundaryEta.clear();
579  stsCP_track_boundaryPhi.clear();
580  stsCP_track_boundaryPx.clear();
581  stsCP_track_boundaryPy.clear();
582  stsCP_track_boundaryPz.clear();
594 
604  simTICLCandidate_time.clear();
606  simTICLCandidate_pdgId.clear();
607  simTICLCandidate_charge.clear();
609 
610  nCandidates = 0;
611  candidate_charge.clear();
612  candidate_pdgId.clear();
613  candidate_energy.clear();
614  candidate_raw_energy.clear();
615  candidate_px.clear();
616  candidate_py.clear();
617  candidate_pz.clear();
618  candidate_time.clear();
619  candidate_time_err.clear();
621  tracksters_in_candidate.clear();
622  track_in_candidate.clear();
623 
624  nTrackstersMerged = 0;
625  tracksters_merged_time.clear();
630  tracksters_merged_raw_pt.clear();
637  tracksters_merged_EV1.clear();
638  tracksters_merged_EV2.clear();
639  tracksters_merged_EV3.clear();
647  tracksters_merged_time.clear();
652  tracksters_merged_raw_pt.clear();
654 
665 
672 
679 
686 
693 
700 
701  nsimTrackstersSC = 0;
702 
703  cluster_seedID.clear();
704  cluster_energy.clear();
705  cluster_correctedEnergy.clear();
707  cluster_position_x.clear();
708  cluster_position_y.clear();
709  cluster_position_z.clear();
710  cluster_position_eta.clear();
711  cluster_position_phi.clear();
712  cluster_layer_id.clear();
713  cluster_type.clear();
714  cluster_time.clear();
715  cluster_timeErr.clear();
716  cluster_number_of_hits.clear();
717 
718  track_id.clear();
719  track_hgcal_x.clear();
720  track_hgcal_y.clear();
721  track_hgcal_z.clear();
722  track_hgcal_eta.clear();
723  track_hgcal_phi.clear();
724  track_hgcal_px.clear();
725  track_hgcal_py.clear();
726  track_hgcal_pz.clear();
727  track_hgcal_pt.clear();
728  track_quality.clear();
729  track_pt.clear();
730  track_missing_outer_hits.clear();
731  track_missing_inner_hits.clear();
732  track_charge.clear();
733  track_time.clear();
734  track_time_quality.clear();
735  track_time_err.clear();
736  track_beta.clear();
737  track_time_mtd.clear();
738  track_time_mtd_err.clear();
739  track_pos_mtd.clear();
740  track_nhits.clear();
741  track_isMuon.clear();
742  track_isTrackerMuon.clear();
743 };
744 
746  : tracksters_token_(consumes<std::vector<ticl::Trackster>>(ps.getParameter<edm::InputTag>("trackstersclue3d"))),
747  tracksters_in_candidate_token_(
748  consumes<std::vector<ticl::Trackster>>(ps.getParameter<edm::InputTag>("trackstersInCand"))),
749  layer_clusters_token_(consumes<std::vector<reco::CaloCluster>>(ps.getParameter<edm::InputTag>("layerClusters"))),
750  ticl_candidates_token_(consumes<std::vector<TICLCandidate>>(ps.getParameter<edm::InputTag>("ticlcandidates"))),
751  tracks_token_(consumes<std::vector<reco::Track>>(ps.getParameter<edm::InputTag>("tracks"))),
752  tracks_time_token_(consumes<edm::ValueMap<float>>(ps.getParameter<edm::InputTag>("tracksTime"))),
753  tracks_time_quality_token_(consumes<edm::ValueMap<float>>(ps.getParameter<edm::InputTag>("tracksTimeQual"))),
754  tracks_time_err_token_(consumes<edm::ValueMap<float>>(ps.getParameter<edm::InputTag>("tracksTimeErr"))),
755  tracks_beta_token_(consumes<edm::ValueMap<float>>(ps.getParameter<edm::InputTag>("tracksBeta"))),
756  tracks_time_mtd_token_(consumes<edm::ValueMap<float>>(ps.getParameter<edm::InputTag>("tracksTimeMtd"))),
757  tracks_time_mtd_err_token_(consumes<edm::ValueMap<float>>(ps.getParameter<edm::InputTag>("tracksTimeMtdErr"))),
758  tracks_pos_mtd_token_(consumes<edm::ValueMap<GlobalPoint>>(ps.getParameter<edm::InputTag>("tracksPosMtd"))),
759  tracksters_merged_token_(
760  consumes<std::vector<ticl::Trackster>>(ps.getParameter<edm::InputTag>("trackstersmerged"))),
761  muons_token_(consumes<std::vector<reco::Muon>>(ps.getParameter<edm::InputTag>("muons"))),
762  clustersTime_token_(
763  consumes<edm::ValueMap<std::pair<float, float>>>(ps.getParameter<edm::InputTag>("layer_clustersTime"))),
764  caloGeometry_token_(esConsumes<CaloGeometry, CaloGeometryRecord, edm::Transition::BeginRun>()),
765  simTracksters_SC_token_(
766  consumes<std::vector<ticl::Trackster>>(ps.getParameter<edm::InputTag>("simtrackstersSC"))),
767  simTracksters_CP_token_(
768  consumes<std::vector<ticl::Trackster>>(ps.getParameter<edm::InputTag>("simtrackstersCP"))),
769  simTracksters_PU_token_(
770  consumes<std::vector<ticl::Trackster>>(ps.getParameter<edm::InputTag>("simtrackstersPU"))),
771  simTICLCandidate_token_(
772  consumes<std::vector<TICLCandidate>>(ps.getParameter<edm::InputTag>("simTICLCandidates"))),
773  tsRecoToSimSC_token_(
774  consumes<ticl::RecoToSimCollectionSimTracksters>(ps.getParameter<edm::InputTag>("recoToSimAssociatorSC"))),
775  tsSimToRecoSC_token_(
776  consumes<ticl::SimToRecoCollectionSimTracksters>(ps.getParameter<edm::InputTag>("simToRecoAssociatorSC"))),
777  tsRecoToSimCP_token_(
778  consumes<ticl::RecoToSimCollectionSimTracksters>(ps.getParameter<edm::InputTag>("recoToSimAssociatorCP"))),
779  tsSimToRecoCP_token_(
780  consumes<ticl::SimToRecoCollectionSimTracksters>(ps.getParameter<edm::InputTag>("simToRecoAssociatorCP"))),
781  MergeRecoToSimSC_token_(consumes<ticl::RecoToSimCollectionSimTracksters>(
782  ps.getParameter<edm::InputTag>("MergerecoToSimAssociatorSC"))),
783  MergeSimToRecoSC_token_(consumes<ticl::SimToRecoCollectionSimTracksters>(
784  ps.getParameter<edm::InputTag>("MergesimToRecoAssociatorSC"))),
785  MergeRecoToSimCP_token_(consumes<ticl::RecoToSimCollectionSimTracksters>(
786  ps.getParameter<edm::InputTag>("MergerecoToSimAssociatorCP"))),
787  MergeSimToRecoCP_token_(consumes<ticl::SimToRecoCollectionSimTracksters>(
788  ps.getParameter<edm::InputTag>("MergesimToRecoAssociatorCP"))),
789  MergeRecoToSimPU_token_(consumes<ticl::RecoToSimCollectionSimTracksters>(
790  ps.getParameter<edm::InputTag>("MergerecoToSimAssociatorPU"))),
791  MergeSimToRecoPU_token_(consumes<ticl::SimToRecoCollectionSimTracksters>(
792  ps.getParameter<edm::InputTag>("MergesimToRecoAssociatorPU"))),
793  simclusters_token_(consumes(ps.getParameter<edm::InputTag>("simclusters"))),
794  caloparticles_token_(consumes(ps.getParameter<edm::InputTag>("caloparticles"))),
795  geometry_token_(esConsumes<CaloGeometry, CaloGeometryRecord, edm::Transition::BeginRun>()),
796  detector_(ps.getParameter<std::string>("detector")),
797  propName_(ps.getParameter<std::string>("propagator")),
798  bfield_token_(esConsumes<MagneticField, IdealMagneticFieldRecord, edm::Transition::BeginRun>()),
799  propagator_token_(
801  saveLCs_(ps.getParameter<bool>("saveLCs")),
802  saveCLUE3DTracksters_(ps.getParameter<bool>("saveCLUE3DTracksters")),
803  saveTrackstersMerged_(ps.getParameter<bool>("saveTrackstersMerged")),
804  saveSimTrackstersSC_(ps.getParameter<bool>("saveSimTrackstersSC")),
805  saveSimTrackstersCP_(ps.getParameter<bool>("saveSimTrackstersCP")),
806  saveTICLCandidate_(ps.getParameter<bool>("saveSimTICLCandidate")),
807  saveSimTICLCandidate_(ps.getParameter<bool>("saveSimTICLCandidate")),
808  saveTracks_(ps.getParameter<bool>("saveTracks")),
809  saveAssociations_(ps.getParameter<bool>("saveAssociations")) {
810  std::string detectorName_ = (detector_ == "HFNose") ? "HGCalHFNoseSensitive" : "HGCalEESensitive";
811  hdc_token_ =
812  esConsumes<HGCalDDDConstants, IdealGeometryRecord, edm::Transition::BeginRun>(edm::ESInputTag("", detectorName_));
813 };
814 
816 
820 
822  hgcons_ = hdc.product();
826 }
827 
828 // Define tree and branches
831  if (saveCLUE3DTracksters_) {
832  trackster_tree_ = fs->make<TTree>("tracksters", "TICL tracksters");
833  trackster_tree_->Branch("event", &ev_event_);
834  trackster_tree_->Branch("NClusters", &nclusters_);
835  trackster_tree_->Branch("NTracksters", &ntracksters_);
836  trackster_tree_->Branch("time", &trackster_time);
837  trackster_tree_->Branch("timeError", &trackster_timeError);
838  trackster_tree_->Branch("regressed_energy", &trackster_regressed_energy);
839  trackster_tree_->Branch("raw_energy", &trackster_raw_energy);
840  trackster_tree_->Branch("raw_em_energy", &trackster_raw_em_energy);
841  trackster_tree_->Branch("raw_pt", &trackster_raw_pt);
842  trackster_tree_->Branch("raw_em_pt", &trackster_raw_em_pt);
843  trackster_tree_->Branch("barycenter_x", &trackster_barycenter_x);
844  trackster_tree_->Branch("barycenter_y", &trackster_barycenter_y);
845  trackster_tree_->Branch("barycenter_z", &trackster_barycenter_z);
846  trackster_tree_->Branch("barycenter_eta", &trackster_barycenter_eta);
847  trackster_tree_->Branch("barycenter_phi", &trackster_barycenter_phi);
848  trackster_tree_->Branch("EV1", &trackster_EV1);
849  trackster_tree_->Branch("EV2", &trackster_EV2);
850  trackster_tree_->Branch("EV3", &trackster_EV3);
851  trackster_tree_->Branch("eVector0_x", &trackster_eVector0_x);
852  trackster_tree_->Branch("eVector0_y", &trackster_eVector0_y);
853  trackster_tree_->Branch("eVector0_z", &trackster_eVector0_z);
854  trackster_tree_->Branch("sigmaPCA1", &trackster_sigmaPCA1);
855  trackster_tree_->Branch("sigmaPCA2", &trackster_sigmaPCA2);
856  trackster_tree_->Branch("sigmaPCA3", &trackster_sigmaPCA3);
857  trackster_tree_->Branch("id_probabilities", &trackster_id_probabilities);
858  trackster_tree_->Branch("vertices_indexes", &trackster_vertices_indexes);
859  trackster_tree_->Branch("vertices_x", &trackster_vertices_x);
860  trackster_tree_->Branch("vertices_y", &trackster_vertices_y);
861  trackster_tree_->Branch("vertices_z", &trackster_vertices_z);
862  trackster_tree_->Branch("vertices_time", &trackster_vertices_time);
863  trackster_tree_->Branch("vertices_timeErr", &trackster_vertices_timeErr);
864  trackster_tree_->Branch("vertices_energy", &trackster_vertices_energy);
865  trackster_tree_->Branch("vertices_correctedEnergy", &trackster_vertices_correctedEnergy);
866  trackster_tree_->Branch("vertices_correctedEnergyUncertainty", &trackster_vertices_correctedEnergyUncertainty);
867  trackster_tree_->Branch("vertices_multiplicity", &trackster_vertices_multiplicity);
868  }
869  if (saveLCs_) {
870  cluster_tree_ = fs->make<TTree>("clusters", "TICL tracksters");
871  cluster_tree_->Branch("seedID", &cluster_seedID);
872  cluster_tree_->Branch("energy", &cluster_energy);
873  cluster_tree_->Branch("correctedEnergy", &cluster_correctedEnergy);
874  cluster_tree_->Branch("correctedEnergyUncertainty", &cluster_correctedEnergyUncertainty);
875  cluster_tree_->Branch("position_x", &cluster_position_x);
876  cluster_tree_->Branch("position_y", &cluster_position_y);
877  cluster_tree_->Branch("position_z", &cluster_position_z);
878  cluster_tree_->Branch("position_eta", &cluster_position_eta);
879  cluster_tree_->Branch("position_phi", &cluster_position_phi);
880  cluster_tree_->Branch("cluster_layer_id", &cluster_layer_id);
881  cluster_tree_->Branch("cluster_type", &cluster_type);
882  cluster_tree_->Branch("cluster_time", &cluster_time);
883  cluster_tree_->Branch("cluster_timeErr", &cluster_timeErr);
884  cluster_tree_->Branch("cluster_number_of_hits", &cluster_number_of_hits);
885  }
886  if (saveTICLCandidate_) {
887  candidate_tree_ = fs->make<TTree>("candidates", "TICL candidates");
888  candidate_tree_->Branch("NCandidates", &nCandidates);
889  candidate_tree_->Branch("candidate_charge", &candidate_charge);
890  candidate_tree_->Branch("candidate_pdgId", &candidate_pdgId);
891  candidate_tree_->Branch("candidate_id_probabilities", &candidate_id_probabilities);
892  candidate_tree_->Branch("candidate_time", &candidate_time);
893  candidate_tree_->Branch("candidate_timeErr", &candidate_time_err);
894  candidate_tree_->Branch("candidate_energy", &candidate_energy);
895  candidate_tree_->Branch("candidate_raw_energy", &candidate_raw_energy);
896  candidate_tree_->Branch("candidate_px", &candidate_px);
897  candidate_tree_->Branch("candidate_py", &candidate_py);
898  candidate_tree_->Branch("candidate_pz", &candidate_pz);
899  candidate_tree_->Branch("track_in_candidate", &track_in_candidate);
900  candidate_tree_->Branch("tracksters_in_candidate", &tracksters_in_candidate);
901  }
902  if (saveTrackstersMerged_) {
903  tracksters_merged_tree_ = fs->make<TTree>("trackstersMerged", "TICL tracksters merged");
904  tracksters_merged_tree_->Branch("event", &ev_event_);
907  tracksters_merged_tree_->Branch("regressed_energy", &tracksters_merged_regressed_energy);
909  tracksters_merged_tree_->Branch("raw_em_energy", &tracksters_merged_raw_em_energy);
912  tracksters_merged_tree_->Branch("NTrackstersMerged", &nTrackstersMerged);
916  tracksters_merged_tree_->Branch("barycenter_eta", &tracksters_merged_barycenter_eta);
917  tracksters_merged_tree_->Branch("barycenter_phi", &tracksters_merged_barycenter_phi);
927  tracksters_merged_tree_->Branch("id_probabilities", &tracksters_merged_id_probabilities);
928  tracksters_merged_tree_->Branch("vertices_indexes", &tracksters_merged_vertices_indexes);
932  tracksters_merged_tree_->Branch("vertices_time", &tracksters_merged_vertices_time);
933  tracksters_merged_tree_->Branch("vertices_timeErr", &tracksters_merged_vertices_timeErr);
934  tracksters_merged_tree_->Branch("vertices_energy", &tracksters_merged_vertices_energy);
935  tracksters_merged_tree_->Branch("vertices_correctedEnergy", &tracksters_merged_vertices_correctedEnergy);
936  tracksters_merged_tree_->Branch("vertices_correctedEnergyUncertainty",
938  tracksters_merged_tree_->Branch("vertices_multiplicity", &tracksters_merged_vertices_multiplicity);
939  }
940  if (saveAssociations_) {
941  associations_tree_ = fs->make<TTree>("associations", "Associations");
942  associations_tree_->Branch("tsCLUE3D_recoToSim_SC", &trackstersCLUE3D_recoToSim_SC);
943  associations_tree_->Branch("tsCLUE3D_recoToSim_SC_score", &trackstersCLUE3D_recoToSim_SC_score);
944  associations_tree_->Branch("tsCLUE3D_recoToSim_SC_sharedE", &trackstersCLUE3D_recoToSim_SC_sharedE);
945  associations_tree_->Branch("tsCLUE3D_simToReco_SC", &trackstersCLUE3D_simToReco_SC);
946  associations_tree_->Branch("tsCLUE3D_simToReco_SC_score", &trackstersCLUE3D_simToReco_SC_score);
947  associations_tree_->Branch("tsCLUE3D_simToReco_SC_sharedE", &trackstersCLUE3D_simToReco_SC_sharedE);
948 
949  associations_tree_->Branch("tsCLUE3D_recoToSim_CP", &trackstersCLUE3D_recoToSim_CP);
950  associations_tree_->Branch("tsCLUE3D_recoToSim_CP_score", &trackstersCLUE3D_recoToSim_CP_score);
951  associations_tree_->Branch("tsCLUE3D_recoToSim_CP_sharedE", &trackstersCLUE3D_recoToSim_CP_sharedE);
952  associations_tree_->Branch("tsCLUE3D_simToReco_CP", &trackstersCLUE3D_simToReco_CP);
953  associations_tree_->Branch("tsCLUE3D_simToReco_CP_score", &trackstersCLUE3D_simToReco_CP_score);
954  associations_tree_->Branch("tsCLUE3D_simToReco_CP_sharedE", &trackstersCLUE3D_simToReco_CP_sharedE);
955 
956  associations_tree_->Branch("Mergetstracksters_recoToSim_SC", &MergeTracksters_recoToSim_SC);
957  associations_tree_->Branch("Mergetstracksters_recoToSim_SC_score", &MergeTracksters_recoToSim_SC_score);
958  associations_tree_->Branch("Mergetstracksters_recoToSim_SC_sharedE", &MergeTracksters_recoToSim_SC_sharedE);
959  associations_tree_->Branch("Mergetstracksters_simToReco_SC", &MergeTracksters_simToReco_SC);
960  associations_tree_->Branch("Mergetstracksters_simToReco_SC_score", &MergeTracksters_simToReco_SC_score);
961  associations_tree_->Branch("Mergetstracksters_simToReco_SC_sharedE", &MergeTracksters_simToReco_SC_sharedE);
962 
963  associations_tree_->Branch("Mergetracksters_recoToSim_CP", &MergeTracksters_recoToSim_CP);
964  associations_tree_->Branch("Mergetracksters_recoToSim_CP_score", &MergeTracksters_recoToSim_CP_score);
965  associations_tree_->Branch("Mergetracksters_recoToSim_CP_sharedE", &MergeTracksters_recoToSim_CP_sharedE);
966  associations_tree_->Branch("Mergetracksters_simToReco_CP", &MergeTracksters_simToReco_CP);
967  associations_tree_->Branch("Mergetracksters_simToReco_CP_score", &MergeTracksters_simToReco_CP_score);
968  associations_tree_->Branch("Mergetracksters_simToReco_CP_sharedE", &MergeTracksters_simToReco_CP_sharedE);
969 
970  associations_tree_->Branch("Mergetracksters_recoToSim_PU", &MergeTracksters_recoToSim_PU);
971  associations_tree_->Branch("Mergetracksters_recoToSim_PU_score", &MergeTracksters_recoToSim_PU_score);
972  associations_tree_->Branch("Mergetracksters_recoToSim_PU_sharedE", &MergeTracksters_recoToSim_PU_sharedE);
973  associations_tree_->Branch("Mergetracksters_simToReco_PU", &MergeTracksters_simToReco_PU);
974  associations_tree_->Branch("Mergetracksters_simToReco_PU_score", &MergeTracksters_simToReco_PU_score);
975  associations_tree_->Branch("Mergetracksters_simToReco_PU_sharedE", &MergeTracksters_simToReco_PU_sharedE);
976  }
977 
978  if (saveSimTrackstersSC_) {
979  simtrackstersSC_tree_ = fs->make<TTree>("simtrackstersSC", "TICL simTracksters SC");
980  simtrackstersSC_tree_->Branch("event", &ev_event_);
981  simtrackstersSC_tree_->Branch("NTracksters", &stsSC_ntracksters_);
983  simtrackstersSC_tree_->Branch("timeBoundary", &stsSC_trackster_timeBoundary);
984  simtrackstersSC_tree_->Branch("timeError", &stsSC_trackster_timeError);
985  simtrackstersSC_tree_->Branch("regressed_energy", &stsSC_trackster_regressed_energy);
986  simtrackstersSC_tree_->Branch("regressed_pt", &stsSC_trackster_regressed_pt);
987  simtrackstersSC_tree_->Branch("raw_energy", &stsSC_trackster_raw_energy);
988  simtrackstersSC_tree_->Branch("raw_em_energy", &stsSC_trackster_raw_em_energy);
989  simtrackstersSC_tree_->Branch("raw_pt", &stsSC_trackster_raw_pt);
990  simtrackstersSC_tree_->Branch("raw_em_pt", &stsSC_trackster_raw_em_pt);
991  simtrackstersSC_tree_->Branch("barycenter_x", &stsSC_trackster_barycenter_x);
992  simtrackstersSC_tree_->Branch("barycenter_y", &stsSC_trackster_barycenter_y);
993  simtrackstersSC_tree_->Branch("barycenter_z", &stsSC_trackster_barycenter_z);
994  simtrackstersSC_tree_->Branch("barycenter_eta", &stsSC_trackster_barycenter_eta);
995  simtrackstersSC_tree_->Branch("barycenter_phi", &stsSC_trackster_barycenter_phi);
999  simtrackstersSC_tree_->Branch("eVector0_x", &stsSC_trackster_eVector0_x);
1000  simtrackstersSC_tree_->Branch("eVector0_y", &stsSC_trackster_eVector0_y);
1001  simtrackstersSC_tree_->Branch("eVector0_z", &stsSC_trackster_eVector0_z);
1002  simtrackstersSC_tree_->Branch("sigmaPCA1", &stsSC_trackster_sigmaPCA1);
1003  simtrackstersSC_tree_->Branch("sigmaPCA2", &stsSC_trackster_sigmaPCA2);
1004  simtrackstersSC_tree_->Branch("sigmaPCA3", &stsSC_trackster_sigmaPCA3);
1005  simtrackstersSC_tree_->Branch("pdgID", &stsSC_pdgID);
1006  simtrackstersSC_tree_->Branch("trackIdx", &stsSC_trackIdx);
1007  simtrackstersSC_tree_->Branch("trackTime", &stsSC_trackTime);
1008  simtrackstersSC_tree_->Branch("boundaryX", &stsSC_boundaryX);
1009  simtrackstersSC_tree_->Branch("boundaryY", &stsSC_boundaryY);
1010  simtrackstersSC_tree_->Branch("boundaryZ", &stsSC_boundaryZ);
1011  simtrackstersSC_tree_->Branch("boundaryEta", &stsSC_boundaryEta);
1012  simtrackstersSC_tree_->Branch("boundaryPhi", &stsSC_boundaryPhi);
1013  simtrackstersSC_tree_->Branch("boundaryPx", &stsSC_boundaryPx);
1014  simtrackstersSC_tree_->Branch("boundaryPy", &stsSC_boundaryPy);
1015  simtrackstersSC_tree_->Branch("boundaryPz", &stsSC_boundaryPz);
1016  simtrackstersSC_tree_->Branch("track_boundaryX", &stsSC_track_boundaryX);
1017  simtrackstersSC_tree_->Branch("track_boundaryY", &stsSC_track_boundaryY);
1018  simtrackstersSC_tree_->Branch("track_boundaryZ", &stsSC_track_boundaryZ);
1019  simtrackstersSC_tree_->Branch("track_boundaryEta", &stsSC_track_boundaryEta);
1020  simtrackstersSC_tree_->Branch("track_boundaryPhi", &stsSC_track_boundaryPhi);
1021  simtrackstersSC_tree_->Branch("track_boundaryPx", &stsSC_track_boundaryPx);
1022  simtrackstersSC_tree_->Branch("track_boundaryPy", &stsSC_track_boundaryPy);
1023  simtrackstersSC_tree_->Branch("track_boundaryPz", &stsSC_track_boundaryPz);
1024  simtrackstersSC_tree_->Branch("id_probabilities", &stsSC_trackster_id_probabilities);
1025  simtrackstersSC_tree_->Branch("vertices_indexes", &stsSC_trackster_vertices_indexes);
1026  simtrackstersSC_tree_->Branch("vertices_x", &stsSC_trackster_vertices_x);
1027  simtrackstersSC_tree_->Branch("vertices_y", &stsSC_trackster_vertices_y);
1028  simtrackstersSC_tree_->Branch("vertices_z", &stsSC_trackster_vertices_z);
1029  simtrackstersSC_tree_->Branch("vertices_time", &stsSC_trackster_vertices_time);
1030  simtrackstersSC_tree_->Branch("vertices_timeErr", &stsSC_trackster_vertices_timeErr);
1031  simtrackstersSC_tree_->Branch("vertices_energy", &stsSC_trackster_vertices_energy);
1032  simtrackstersSC_tree_->Branch("vertices_correctedEnergy", &stsSC_trackster_vertices_correctedEnergy);
1033  simtrackstersSC_tree_->Branch("vertices_correctedEnergyUncertainty",
1035  simtrackstersSC_tree_->Branch("vertices_multiplicity", &stsSC_trackster_vertices_multiplicity);
1036  simtrackstersSC_tree_->Branch("NsimTrackstersSC", &nsimTrackstersSC);
1037  }
1038  if (saveSimTrackstersCP_) {
1039  simtrackstersCP_tree_ = fs->make<TTree>("simtrackstersCP", "TICL simTracksters CP");
1040  simtrackstersCP_tree_->Branch("event", &ev_event_);
1041  simtrackstersCP_tree_->Branch("NTracksters", &stsCP_ntracksters_);
1042  simtrackstersCP_tree_->Branch("time", &stsCP_trackster_time);
1043  simtrackstersCP_tree_->Branch("timeBoundary", &stsCP_trackster_timeBoundary);
1044  simtrackstersCP_tree_->Branch("timeError", &stsCP_trackster_timeError);
1045  simtrackstersCP_tree_->Branch("regressed_energy", &stsCP_trackster_regressed_energy);
1046  simtrackstersCP_tree_->Branch("regressed_pt", &stsCP_trackster_regressed_pt);
1047  simtrackstersCP_tree_->Branch("raw_energy", &stsCP_trackster_raw_energy);
1048  simtrackstersCP_tree_->Branch("raw_em_energy", &stsCP_trackster_raw_em_energy);
1049  simtrackstersCP_tree_->Branch("raw_pt", &stsCP_trackster_raw_pt);
1050  simtrackstersCP_tree_->Branch("raw_em_pt", &stsCP_trackster_raw_em_pt);
1051  simtrackstersCP_tree_->Branch("barycenter_x", &stsCP_trackster_barycenter_x);
1052  simtrackstersCP_tree_->Branch("barycenter_y", &stsCP_trackster_barycenter_y);
1053  simtrackstersCP_tree_->Branch("barycenter_z", &stsCP_trackster_barycenter_z);
1054  simtrackstersCP_tree_->Branch("barycenter_eta", &stsCP_trackster_barycenter_eta);
1055  simtrackstersCP_tree_->Branch("barycenter_phi", &stsCP_trackster_barycenter_phi);
1056  simtrackstersCP_tree_->Branch("pdgID", &stsCP_pdgID);
1057  simtrackstersCP_tree_->Branch("trackIdx", &stsCP_trackIdx);
1058  simtrackstersCP_tree_->Branch("trackTime", &stsCP_trackTime);
1059  simtrackstersCP_tree_->Branch("boundaryX", &stsCP_boundaryX);
1060  simtrackstersCP_tree_->Branch("boundaryY", &stsCP_boundaryY);
1061  simtrackstersCP_tree_->Branch("boundaryZ", &stsCP_boundaryZ);
1062  simtrackstersCP_tree_->Branch("boundaryEta", &stsCP_boundaryEta);
1063  simtrackstersCP_tree_->Branch("boundaryPhi", &stsCP_boundaryPhi);
1064  simtrackstersCP_tree_->Branch("boundaryPx", &stsCP_boundaryPx);
1065  simtrackstersCP_tree_->Branch("boundaryPy", &stsCP_boundaryPy);
1066  simtrackstersCP_tree_->Branch("boundaryPz", &stsCP_boundaryPz);
1067  simtrackstersCP_tree_->Branch("track_boundaryX", &stsCP_track_boundaryX);
1068  simtrackstersCP_tree_->Branch("track_boundaryY", &stsCP_track_boundaryY);
1069  simtrackstersCP_tree_->Branch("track_boundaryZ", &stsCP_track_boundaryZ);
1070  simtrackstersCP_tree_->Branch("track_boundaryEta", &stsCP_track_boundaryEta);
1071  simtrackstersCP_tree_->Branch("track_boundaryPhi", &stsCP_track_boundaryPhi);
1072  simtrackstersCP_tree_->Branch("track_boundaryPx", &stsCP_track_boundaryPx);
1073  simtrackstersCP_tree_->Branch("track_boundaryPy", &stsCP_track_boundaryPy);
1074  simtrackstersCP_tree_->Branch("track_boundaryPz", &stsCP_track_boundaryPz);
1075  simtrackstersCP_tree_->Branch("EV1", &stsCP_trackster_EV1);
1076  simtrackstersCP_tree_->Branch("EV2", &stsCP_trackster_EV2);
1077  simtrackstersCP_tree_->Branch("EV3", &stsCP_trackster_EV3);
1078  simtrackstersCP_tree_->Branch("eVector0_x", &stsCP_trackster_eVector0_x);
1079  simtrackstersCP_tree_->Branch("eVector0_y", &stsCP_trackster_eVector0_y);
1080  simtrackstersCP_tree_->Branch("eVector0_z", &stsCP_trackster_eVector0_z);
1081  simtrackstersCP_tree_->Branch("sigmaPCA1", &stsCP_trackster_sigmaPCA1);
1082  simtrackstersCP_tree_->Branch("sigmaPCA2", &stsCP_trackster_sigmaPCA2);
1083  simtrackstersCP_tree_->Branch("sigmaPCA3", &stsCP_trackster_sigmaPCA3);
1084  simtrackstersCP_tree_->Branch("id_probabilities", &stsCP_trackster_id_probabilities);
1085  simtrackstersCP_tree_->Branch("vertices_indexes", &stsCP_trackster_vertices_indexes);
1086  simtrackstersCP_tree_->Branch("vertices_x", &stsCP_trackster_vertices_x);
1087  simtrackstersCP_tree_->Branch("vertices_y", &stsCP_trackster_vertices_y);
1088  simtrackstersCP_tree_->Branch("vertices_z", &stsCP_trackster_vertices_z);
1089  simtrackstersCP_tree_->Branch("vertices_time", &stsCP_trackster_vertices_time);
1090  simtrackstersCP_tree_->Branch("vertices_timeErr", &stsCP_trackster_vertices_timeErr);
1091  simtrackstersCP_tree_->Branch("vertices_energy", &stsCP_trackster_vertices_energy);
1092  simtrackstersCP_tree_->Branch("vertices_correctedEnergy", &stsCP_trackster_vertices_correctedEnergy);
1093  simtrackstersCP_tree_->Branch("vertices_correctedEnergyUncertainty",
1095  simtrackstersCP_tree_->Branch("vertices_multiplicity", &stsCP_trackster_vertices_multiplicity);
1096  }
1097 
1098  if (saveTracks_) {
1099  tracks_tree_ = fs->make<TTree>("tracks", "Tracks");
1100  tracks_tree_->Branch("event", &ev_event_);
1101  tracks_tree_->Branch("track_id", &track_id);
1102  tracks_tree_->Branch("track_hgcal_x", &track_hgcal_x);
1103  tracks_tree_->Branch("track_hgcal_y", &track_hgcal_y);
1104  tracks_tree_->Branch("track_hgcal_z", &track_hgcal_z);
1105  tracks_tree_->Branch("track_hgcal_eta", &track_hgcal_eta);
1106  tracks_tree_->Branch("track_hgcal_phi", &track_hgcal_phi);
1107  tracks_tree_->Branch("track_hgcal_pt", &track_hgcal_pt);
1108  tracks_tree_->Branch("track_pt", &track_pt);
1109  tracks_tree_->Branch("track_missing_outer_hits", &track_missing_outer_hits);
1110  tracks_tree_->Branch("track_missing_inner_hits", &track_missing_inner_hits);
1111  tracks_tree_->Branch("track_quality", &track_quality);
1112  tracks_tree_->Branch("track_charge", &track_charge);
1113  tracks_tree_->Branch("track_time", &track_time);
1114  tracks_tree_->Branch("track_time_quality", &track_time_quality);
1115  tracks_tree_->Branch("track_time_err", &track_time_err);
1116  tracks_tree_->Branch("track_beta", &track_beta);
1117  tracks_tree_->Branch("track_time_mtd", &track_time_mtd);
1118  tracks_tree_->Branch("track_time_mtd_err", &track_time_mtd_err);
1119  tracks_tree_->Branch("track_pos_mtd", &track_pos_mtd);
1120  tracks_tree_->Branch("track_nhits", &track_nhits);
1121  tracks_tree_->Branch("track_isMuon", &track_isMuon);
1122  tracks_tree_->Branch("track_isTrackerMuon", &track_isTrackerMuon);
1123  }
1124 
1125  if (saveSimTICLCandidate_) {
1126  simTICLCandidate_tree = fs->make<TTree>("simTICLCandidate", "Sim TICL Candidate");
1127  simTICLCandidate_tree->Branch("simTICLCandidate_raw_energy", &simTICLCandidate_raw_energy);
1128  simTICLCandidate_tree->Branch("simTICLCandidate_regressed_energy", &simTICLCandidate_regressed_energy);
1129  simTICLCandidate_tree->Branch("simTICLCandidate_simTracksterCPIndex", &simTICLCandidate_simTracksterCPIndex);
1130  simTICLCandidate_tree->Branch("simTICLCandidate_boundaryX", &simTICLCandidate_boundaryX);
1131  simTICLCandidate_tree->Branch("simTICLCandidate_boundaryY", &simTICLCandidate_boundaryY);
1132  simTICLCandidate_tree->Branch("simTICLCandidate_boundaryZ", &simTICLCandidate_boundaryZ);
1133  simTICLCandidate_tree->Branch("simTICLCandidate_boundaryPx", &simTICLCandidate_boundaryPx);
1134  simTICLCandidate_tree->Branch("simTICLCandidate_boundaryPy", &simTICLCandidate_boundaryPy);
1135  simTICLCandidate_tree->Branch("simTICLCandidate_boundaryPz", &simTICLCandidate_boundaryPz);
1136  simTICLCandidate_tree->Branch("simTICLCandidate_time", &simTICLCandidate_time);
1137  simTICLCandidate_tree->Branch("simTICLCandidate_caloParticleMass", &simTICLCandidate_caloParticleMass);
1138  simTICLCandidate_tree->Branch("simTICLCandidate_pdgId", &simTICLCandidate_pdgId);
1139  simTICLCandidate_tree->Branch("simTICLCandidate_charge", &simTICLCandidate_charge);
1140  simTICLCandidate_tree->Branch("simTICLCandidate_track_in_candidate", &simTICLCandidate_track_in_candidate);
1141  }
1142 }
1143 
1145  // build disks at HGCal front & EM-Had interface for track propagation
1146 
1147  float zVal = hgcons_->waferZ(1, true);
1148  std::pair<float, float> rMinMax = hgcons_->rangeR(zVal, true);
1149 
1150  float zVal_interface = rhtools_.getPositionLayer(rhtools_.lastLayerEE()).z();
1151  std::pair<float, float> rMinMax_interface = hgcons_->rangeR(zVal_interface, true);
1152 
1153  for (int iSide = 0; iSide < 2; ++iSide) {
1154  float zSide = (iSide == 0) ? (-1. * zVal) : zVal;
1155  firstDisk_[iSide] =
1156  std::make_unique<GeomDet>(Disk::build(Disk::PositionType(0, 0, zSide),
1158  SimpleDiskBounds(rMinMax.first, rMinMax.second, zSide - 0.5, zSide + 0.5))
1159  .get());
1160 
1161  zSide = (iSide == 0) ? (-1. * zVal_interface) : zVal_interface;
1162  interfaceDisk_[iSide] = std::make_unique<GeomDet>(
1163  Disk::build(Disk::PositionType(0, 0, zSide),
1165  SimpleDiskBounds(rMinMax_interface.first, rMinMax_interface.second, zSide - 0.5, zSide + 0.5))
1166  .get());
1167  }
1168 }
1169 
1171  const hgcal::RecHitTools rhtools,
1172  const edm::ESHandle<MagneticField> bfieldH,
1173  const edm::ESHandle<Propagator> propH) {
1174  hgcons_ = hgcons;
1175  rhtools_ = rhtools;
1176  buildLayers();
1177 
1178  bfield_ = bfieldH;
1179  propagator_ = propH;
1180 }
1182  ev_event_ += 1;
1183  clearVariables();
1184  auto bFieldProd = bfield_.product();
1185  const Propagator& prop = (*propagator_);
1186  //get all the tracksters
1187  edm::Handle<std::vector<ticl::Trackster>> tracksters_handle;
1188  event.getByToken(tracksters_token_, tracksters_handle);
1189  const auto& tracksters = *tracksters_handle;
1190 
1191  edm::Handle<std::vector<ticl::Trackster>> tracksters_in_candidate_handle;
1192  event.getByToken(tracksters_in_candidate_token_, tracksters_in_candidate_handle);
1193 
1194  //get all the layer clusters
1196  event.getByToken(layer_clusters_token_, layer_clusters_h);
1197  const auto& clusters = *layer_clusters_h;
1198 
1200  event.getByToken(clustersTime_token_, clustersTime_h);
1201  const auto& layerClustersTimes = *clustersTime_h;
1202 
1203  //TICL Candidate
1205  event.getByToken(ticl_candidates_token_, candidates_h);
1206  const auto& ticlcandidates = *candidates_h;
1207 
1208  //Track
1210  event.getByToken(tracks_token_, tracks_h);
1211  const auto& tracks = *tracks_h;
1212 
1213  edm::Handle<edm::ValueMap<float>> trackTime_h;
1214  event.getByToken(tracks_time_token_, trackTime_h);
1215  const auto& trackTime = *trackTime_h;
1216 
1217  edm::Handle<edm::ValueMap<float>> trackTimeErr_h;
1218  event.getByToken(tracks_time_err_token_, trackTimeErr_h);
1219  const auto& trackTimeErr = *trackTimeErr_h;
1220 
1221  edm::Handle<edm::ValueMap<float>> trackBeta_h;
1222  event.getByToken(tracks_beta_token_, trackBeta_h);
1223  const auto& trackBeta = *trackBeta_h;
1224 
1225  edm::Handle<edm::ValueMap<float>> trackTimeQual_h;
1226  event.getByToken(tracks_time_quality_token_, trackTimeQual_h);
1227  const auto& trackTimeQual = *trackTimeQual_h;
1228 
1229  edm::Handle<edm::ValueMap<float>> trackTimeMtd_h;
1230  event.getByToken(tracks_time_mtd_token_, trackTimeMtd_h);
1231  const auto& trackTimeMtd = *trackTimeMtd_h;
1232 
1233  edm::Handle<edm::ValueMap<float>> trackTimeMtdErr_h;
1234  event.getByToken(tracks_time_mtd_err_token_, trackTimeMtdErr_h);
1235  const auto& trackTimeMtdErr = *trackTimeMtdErr_h;
1236 
1238  event.getByToken(tracks_pos_mtd_token_, trackPosMtd_h);
1239  const auto& trackPosMtd = *trackPosMtd_h;
1240 
1241  //Tracksters merged
1242  edm::Handle<std::vector<ticl::Trackster>> tracksters_merged_h;
1243  event.getByToken(tracksters_merged_token_, tracksters_merged_h);
1244  const auto& trackstersmerged = *tracksters_merged_h;
1245 
1246  // muons
1248  event.getByToken(muons_token_, muons_h);
1249  auto& muons = *muons_h;
1250 
1251  // simTracksters from SC
1252  edm::Handle<std::vector<ticl::Trackster>> simTrackstersSC_h;
1253  event.getByToken(simTracksters_SC_token_, simTrackstersSC_h);
1254  const auto& simTrackstersSC = *simTrackstersSC_h;
1255 
1256  // simTracksters from CP
1257  edm::Handle<std::vector<ticl::Trackster>> simTrackstersCP_h;
1258  event.getByToken(simTracksters_CP_token_, simTrackstersCP_h);
1259  const auto& simTrackstersCP = *simTrackstersCP_h;
1260 
1261  // simTracksters from PU
1262  edm::Handle<std::vector<ticl::Trackster>> simTrackstersPU_h;
1263  event.getByToken(simTracksters_PU_token_, simTrackstersPU_h);
1264  const auto& simTrackstersPU = *simTrackstersPU_h;
1265 
1266  edm::Handle<std::vector<TICLCandidate>> simTICLCandidates_h;
1267  event.getByToken(simTICLCandidate_token_, simTICLCandidates_h);
1268  const auto& simTICLCandidates = *simTICLCandidates_h;
1269 
1270  // trackster reco to sim SC
1272  event.getByToken(tsRecoToSimSC_token_, tsRecoToSimSC_h);
1273  auto const& tsRecoSimSCMap = *tsRecoToSimSC_h;
1274 
1275  // sim simTrackster SC to reco trackster
1277  event.getByToken(tsSimToRecoSC_token_, tsSimToRecoSC_h);
1278  auto const& tsSimToRecoSCMap = *tsSimToRecoSC_h;
1279 
1280  // trackster reco to sim CP
1282  event.getByToken(tsRecoToSimCP_token_, tsRecoToSimCP_h);
1283  auto const& tsRecoSimCPMap = *tsRecoToSimCP_h;
1284 
1285  // sim simTrackster CP to reco trackster
1287  event.getByToken(tsSimToRecoCP_token_, tsSimToRecoCP_h);
1288  auto const& tsSimToRecoCPMap = *tsSimToRecoCP_h;
1289 
1291  event.getByToken(MergeRecoToSimSC_token_, mergetsRecoToSimSC_h);
1292  auto const& MergetsRecoSimSCMap = *mergetsRecoToSimSC_h;
1293 
1294  // sim simTrackster SC to reco trackster
1296  event.getByToken(MergeSimToRecoSC_token_, mergetsSimToRecoSC_h);
1297  auto const& MergetsSimToRecoSCMap = *mergetsSimToRecoSC_h;
1298 
1299  // trackster reco to sim CP
1301  event.getByToken(MergeRecoToSimCP_token_, mergetsRecoToSimCP_h);
1302  auto const& MergetsRecoSimCPMap = *mergetsRecoToSimCP_h;
1303 
1304  // sim simTrackster CP to reco trackster
1306  event.getByToken(MergeSimToRecoCP_token_, mergetsSimToRecoCP_h);
1307  auto const& MergetsSimToRecoCPMap = *mergetsSimToRecoCP_h;
1308 
1309  // trackster reco to sim PU
1311  event.getByToken(MergeRecoToSimPU_token_, mergetsRecoToSimPU_h);
1312  auto const& MergetsRecoSimPUMap = *mergetsRecoToSimPU_h;
1313 
1314  // sim simTrackster PU to reco trackster
1316  event.getByToken(MergeSimToRecoPU_token_, mergetsSimToRecoPU_h);
1317  auto const& MergetsSimToRecoPUMap = *mergetsSimToRecoPU_h;
1318 
1319  edm::Handle<std::vector<CaloParticle>> caloparticles_h;
1320  event.getByToken(caloparticles_token_, caloparticles_h);
1321  const auto& caloparticles = *caloparticles_h;
1322 
1323  const auto& simclusters = event.get(simclusters_token_);
1324 
1325  ntracksters_ = tracksters.size();
1326  nclusters_ = clusters.size();
1327 
1328  for (auto trackster_iterator = tracksters.begin(); trackster_iterator != tracksters.end(); ++trackster_iterator) {
1329  //per-trackster analysis
1330  trackster_time.push_back(trackster_iterator->time());
1331  trackster_timeError.push_back(trackster_iterator->timeError());
1332  trackster_regressed_energy.push_back(trackster_iterator->regressed_energy());
1333  trackster_raw_energy.push_back(trackster_iterator->raw_energy());
1334  trackster_raw_em_energy.push_back(trackster_iterator->raw_em_energy());
1335  trackster_raw_pt.push_back(trackster_iterator->raw_pt());
1336  trackster_raw_em_pt.push_back(trackster_iterator->raw_em_pt());
1337  trackster_barycenter_x.push_back(trackster_iterator->barycenter().x());
1338  trackster_barycenter_y.push_back(trackster_iterator->barycenter().y());
1339  trackster_barycenter_z.push_back(trackster_iterator->barycenter().z());
1340  trackster_barycenter_eta.push_back(trackster_iterator->barycenter().eta());
1341  trackster_barycenter_phi.push_back(trackster_iterator->barycenter().phi());
1342  trackster_EV1.push_back(trackster_iterator->eigenvalues()[0]);
1343  trackster_EV2.push_back(trackster_iterator->eigenvalues()[1]);
1344  trackster_EV3.push_back(trackster_iterator->eigenvalues()[2]);
1345  trackster_eVector0_x.push_back((trackster_iterator->eigenvectors()[0]).x());
1346  trackster_eVector0_y.push_back((trackster_iterator->eigenvectors()[0]).y());
1347  trackster_eVector0_z.push_back((trackster_iterator->eigenvectors()[0]).z());
1348  trackster_sigmaPCA1.push_back(trackster_iterator->sigmasPCA()[0]);
1349  trackster_sigmaPCA2.push_back(trackster_iterator->sigmasPCA()[1]);
1350  trackster_sigmaPCA3.push_back(trackster_iterator->sigmasPCA()[2]);
1351  std::vector<float> id_probs;
1352  for (size_t i = 0; i < 8; i++)
1353  id_probs.push_back(trackster_iterator->id_probabilities(i));
1354  trackster_id_probabilities.push_back(id_probs);
1355 
1356  // Clusters
1357  std::vector<uint32_t> vertices_indexes;
1358  std::vector<float> vertices_x;
1359  std::vector<float> vertices_y;
1360  std::vector<float> vertices_z;
1361  std::vector<float> vertices_time;
1362  std::vector<float> vertices_timeErr;
1363  std::vector<float> vertices_energy;
1364  std::vector<float> vertices_correctedEnergy;
1365  std::vector<float> vertices_correctedEnergyUncertainty;
1366  for (auto idx : trackster_iterator->vertices()) {
1367  vertices_indexes.push_back(idx);
1368  auto associated_cluster = (*layer_clusters_h)[idx];
1369  vertices_x.push_back(associated_cluster.x());
1370  vertices_y.push_back(associated_cluster.y());
1371  vertices_z.push_back(associated_cluster.z());
1372  vertices_energy.push_back(associated_cluster.energy());
1373  vertices_correctedEnergy.push_back(associated_cluster.correctedEnergy());
1374  vertices_correctedEnergyUncertainty.push_back(associated_cluster.correctedEnergyUncertainty());
1375  vertices_time.push_back(layerClustersTimes.get(idx).first);
1376  vertices_timeErr.push_back(layerClustersTimes.get(idx).second);
1377  }
1378  trackster_vertices_indexes.push_back(vertices_indexes);
1379  trackster_vertices_x.push_back(vertices_x);
1380  trackster_vertices_y.push_back(vertices_y);
1381  trackster_vertices_z.push_back(vertices_z);
1382  trackster_vertices_time.push_back(vertices_time);
1383  trackster_vertices_timeErr.push_back(vertices_timeErr);
1384  trackster_vertices_energy.push_back(vertices_energy);
1385  trackster_vertices_correctedEnergy.push_back(vertices_correctedEnergy);
1386  trackster_vertices_correctedEnergyUncertainty.push_back(vertices_correctedEnergyUncertainty);
1387 
1388  // Multiplicity
1389  std::vector<float> vertices_multiplicity;
1390  for (auto multiplicity : trackster_iterator->vertex_multiplicity()) {
1391  vertices_multiplicity.push_back(multiplicity);
1392  }
1393  trackster_vertices_multiplicity.push_back(vertices_multiplicity);
1394  }
1395 
1396  stsSC_ntracksters_ = simTrackstersSC.size();
1397  using CaloObjectVariant = std::variant<CaloParticle, SimCluster>;
1398  for (auto trackster_iterator = simTrackstersSC.begin(); trackster_iterator != simTrackstersSC.end();
1399  ++trackster_iterator) {
1400  //per-trackster analysis
1401  stsSC_trackster_time.push_back(trackster_iterator->time());
1402  stsSC_trackster_timeBoundary.push_back(trackster_iterator->boundaryTime());
1403  stsSC_trackster_timeError.push_back(trackster_iterator->timeError());
1404  stsSC_trackster_regressed_energy.push_back(trackster_iterator->regressed_energy());
1405  stsSC_trackster_raw_energy.push_back(trackster_iterator->raw_energy());
1406  stsSC_trackster_raw_em_energy.push_back(trackster_iterator->raw_em_energy());
1407  stsSC_trackster_raw_pt.push_back(trackster_iterator->raw_pt());
1408  stsSC_trackster_raw_em_pt.push_back(trackster_iterator->raw_em_pt());
1409  stsSC_trackster_barycenter_x.push_back(trackster_iterator->barycenter().x());
1410  stsSC_trackster_barycenter_y.push_back(trackster_iterator->barycenter().y());
1411  stsSC_trackster_barycenter_z.push_back(trackster_iterator->barycenter().z());
1412  stsSC_trackster_barycenter_eta.push_back(trackster_iterator->barycenter().eta());
1413  stsSC_trackster_barycenter_phi.push_back(trackster_iterator->barycenter().phi());
1414  stsSC_trackster_EV1.push_back(trackster_iterator->eigenvalues()[0]);
1415  stsSC_trackster_EV2.push_back(trackster_iterator->eigenvalues()[1]);
1416  stsSC_trackster_EV3.push_back(trackster_iterator->eigenvalues()[2]);
1417  stsSC_trackster_eVector0_x.push_back((trackster_iterator->eigenvectors()[0]).x());
1418  stsSC_trackster_eVector0_y.push_back((trackster_iterator->eigenvectors()[0]).y());
1419  stsSC_trackster_eVector0_z.push_back((trackster_iterator->eigenvectors()[0]).z());
1420  stsSC_trackster_sigmaPCA1.push_back(trackster_iterator->sigmasPCA()[0]);
1421  stsSC_trackster_sigmaPCA2.push_back(trackster_iterator->sigmasPCA()[1]);
1422  stsSC_trackster_sigmaPCA3.push_back(trackster_iterator->sigmasPCA()[2]);
1423  stsSC_pdgID.push_back(simclusters[trackster_iterator->seedIndex()].pdgId());
1424 
1425  CaloObjectVariant caloObj;
1426  if (trackster_iterator->seedID() == caloparticles_h.id()) {
1427  caloObj = caloparticles[trackster_iterator->seedIndex()];
1428  } else {
1429  caloObj = simclusters[trackster_iterator->seedIndex()];
1430  }
1431 
1432  auto const& simTrack = std::visit([](auto&& obj) { return obj.g4Tracks()[0]; }, caloObj);
1433  auto const& caloPt = std::visit([](auto&& obj) { return obj.pt(); }, caloObj);
1434  stsSC_trackster_regressed_pt.push_back(caloPt);
1435  if (simTrack.crossedBoundary()) {
1436  stsSC_boundaryX.push_back(simTrack.getPositionAtBoundary().x());
1437  stsSC_boundaryY.push_back(simTrack.getPositionAtBoundary().y());
1438  stsSC_boundaryZ.push_back(simTrack.getPositionAtBoundary().z());
1439  stsSC_boundaryEta.push_back(simTrack.getPositionAtBoundary().eta());
1440  stsSC_boundaryPhi.push_back(simTrack.getPositionAtBoundary().phi());
1441  stsSC_boundaryPx.push_back(simTrack.getMomentumAtBoundary().x());
1442  stsSC_boundaryPy.push_back(simTrack.getMomentumAtBoundary().y());
1443  stsSC_boundaryPz.push_back(simTrack.getMomentumAtBoundary().z());
1444  } else {
1445  stsSC_boundaryX.push_back(-999);
1446  stsSC_boundaryY.push_back(-999);
1447  stsSC_boundaryZ.push_back(-999);
1448  stsSC_boundaryEta.push_back(-999);
1449  stsSC_boundaryPhi.push_back(-999);
1450  stsSC_boundaryPx.push_back(-999);
1451  stsSC_boundaryPy.push_back(-999);
1452  stsSC_boundaryPz.push_back(-999);
1453  }
1454  auto const trackIdx = trackster_iterator->trackIdx();
1455  stsSC_trackIdx.push_back(trackIdx);
1456  if (trackIdx != -1) {
1457  const auto& track = tracks[trackIdx];
1458 
1459  int iSide = int(track.eta() > 0);
1460 
1461  const auto& fts = trajectoryStateTransform::outerFreeState((track), bFieldProd);
1462  // to the HGCal front
1463  const auto& tsos = prop.propagate(fts, firstDisk_[iSide]->surface());
1464  if (tsos.isValid()) {
1465  const auto& globalPos = tsos.globalPosition();
1466  const auto& globalMom = tsos.globalMomentum();
1467  stsSC_track_boundaryX.push_back(globalPos.x());
1468  stsSC_track_boundaryY.push_back(globalPos.y());
1469  stsSC_track_boundaryZ.push_back(globalPos.z());
1470  stsSC_track_boundaryEta.push_back(globalPos.eta());
1471  stsSC_track_boundaryPhi.push_back(globalPos.phi());
1472  stsSC_track_boundaryPx.push_back(globalMom.x());
1473  stsSC_track_boundaryPy.push_back(globalMom.y());
1474  stsSC_track_boundaryPz.push_back(globalMom.z());
1475  stsSC_trackTime.push_back(track.t0());
1476  } else {
1477  stsSC_track_boundaryX.push_back(-999);
1478  stsSC_track_boundaryY.push_back(-999);
1479  stsSC_track_boundaryZ.push_back(-999);
1480  stsSC_track_boundaryEta.push_back(-999);
1481  stsSC_track_boundaryPhi.push_back(-999);
1482  stsSC_track_boundaryPx.push_back(-999);
1483  stsSC_track_boundaryPy.push_back(-999);
1484  stsSC_track_boundaryPz.push_back(-999);
1485  }
1486  } else {
1487  stsSC_track_boundaryX.push_back(-999);
1488  stsSC_track_boundaryY.push_back(-999);
1489  stsSC_track_boundaryZ.push_back(-999);
1490  stsSC_track_boundaryEta.push_back(-999);
1491  stsSC_track_boundaryPhi.push_back(-999);
1492  stsSC_track_boundaryPx.push_back(-999);
1493  stsSC_track_boundaryPy.push_back(-999);
1494  stsSC_track_boundaryPz.push_back(-999);
1495  }
1496 
1497  std::vector<float> id_probs;
1498  for (size_t i = 0; i < 8; i++)
1499  id_probs.push_back(trackster_iterator->id_probabilities(i));
1500  stsSC_trackster_id_probabilities.push_back(id_probs);
1501 
1502  // Clusters
1503  std::vector<uint32_t> vertices_indexes;
1504  std::vector<float> vertices_x;
1505  std::vector<float> vertices_y;
1506  std::vector<float> vertices_z;
1507  std::vector<float> vertices_time;
1508  std::vector<float> vertices_timeErr;
1509  std::vector<float> vertices_energy;
1510  std::vector<float> vertices_correctedEnergy;
1511  std::vector<float> vertices_correctedEnergyUncertainty;
1512  for (auto idx : trackster_iterator->vertices()) {
1513  vertices_indexes.push_back(idx);
1514  auto associated_cluster = (*layer_clusters_h)[idx];
1515  vertices_x.push_back(associated_cluster.x());
1516  vertices_y.push_back(associated_cluster.y());
1517  vertices_z.push_back(associated_cluster.z());
1518  vertices_energy.push_back(associated_cluster.energy());
1519  vertices_correctedEnergy.push_back(associated_cluster.correctedEnergy());
1520  vertices_correctedEnergyUncertainty.push_back(associated_cluster.correctedEnergyUncertainty());
1521  vertices_time.push_back(layerClustersTimes.get(idx).first);
1522  vertices_timeErr.push_back(layerClustersTimes.get(idx).second);
1523  }
1524  stsSC_trackster_vertices_indexes.push_back(vertices_indexes);
1525  stsSC_trackster_vertices_x.push_back(vertices_x);
1526  stsSC_trackster_vertices_y.push_back(vertices_y);
1527  stsSC_trackster_vertices_z.push_back(vertices_z);
1528  stsSC_trackster_vertices_time.push_back(vertices_time);
1529  stsSC_trackster_vertices_timeErr.push_back(vertices_timeErr);
1530  stsSC_trackster_vertices_energy.push_back(vertices_energy);
1531  stsSC_trackster_vertices_correctedEnergy.push_back(vertices_correctedEnergy);
1532  stsSC_trackster_vertices_correctedEnergyUncertainty.push_back(vertices_correctedEnergyUncertainty);
1533 
1534  // Multiplicity
1535  std::vector<float> vertices_multiplicity;
1536  for (auto multiplicity : trackster_iterator->vertex_multiplicity()) {
1537  vertices_multiplicity.push_back(multiplicity);
1538  }
1539  stsSC_trackster_vertices_multiplicity.push_back(vertices_multiplicity);
1540  }
1541 
1542  stsCP_ntracksters_ = simTrackstersCP.size();
1543 
1544  for (auto trackster_iterator = simTrackstersCP.begin(); trackster_iterator != simTrackstersCP.end();
1545  ++trackster_iterator) {
1546  //per-trackster analysis
1547  stsCP_trackster_time.push_back(trackster_iterator->time());
1548  stsCP_trackster_timeBoundary.push_back(trackster_iterator->boundaryTime());
1549  stsCP_trackster_timeError.push_back(trackster_iterator->timeError());
1550  stsCP_trackster_regressed_energy.push_back(trackster_iterator->regressed_energy());
1551  stsCP_trackster_raw_energy.push_back(trackster_iterator->raw_energy());
1552  stsCP_trackster_raw_em_energy.push_back(trackster_iterator->raw_em_energy());
1553  stsCP_trackster_raw_pt.push_back(trackster_iterator->raw_pt());
1554  stsCP_trackster_raw_em_pt.push_back(trackster_iterator->raw_em_pt());
1555  stsCP_trackster_barycenter_x.push_back(trackster_iterator->barycenter().x());
1556  stsCP_trackster_barycenter_y.push_back(trackster_iterator->barycenter().y());
1557  stsCP_trackster_barycenter_z.push_back(trackster_iterator->barycenter().z());
1558  stsCP_trackster_barycenter_eta.push_back(trackster_iterator->barycenter().eta());
1559  stsCP_trackster_barycenter_phi.push_back(trackster_iterator->barycenter().phi());
1560  stsCP_trackster_EV1.push_back(trackster_iterator->eigenvalues()[0]);
1561  stsCP_trackster_EV2.push_back(trackster_iterator->eigenvalues()[1]);
1562  stsCP_trackster_EV3.push_back(trackster_iterator->eigenvalues()[2]);
1563  stsCP_trackster_eVector0_x.push_back((trackster_iterator->eigenvectors()[0]).x());
1564  stsCP_trackster_eVector0_y.push_back((trackster_iterator->eigenvectors()[0]).y());
1565  stsCP_trackster_eVector0_z.push_back((trackster_iterator->eigenvectors()[0]).z());
1566  stsCP_trackster_sigmaPCA1.push_back(trackster_iterator->sigmasPCA()[0]);
1567  stsCP_trackster_sigmaPCA2.push_back(trackster_iterator->sigmasPCA()[1]);
1568  stsCP_trackster_sigmaPCA3.push_back(trackster_iterator->sigmasPCA()[2]);
1569  stsCP_pdgID.push_back(caloparticles[trackster_iterator->seedIndex()].pdgId());
1570  CaloObjectVariant caloObj;
1571  if (trackster_iterator->seedID() == caloparticles_h.id()) {
1572  caloObj = caloparticles[trackster_iterator->seedIndex()];
1573  } else {
1574  caloObj = simclusters[trackster_iterator->seedIndex()];
1575  }
1576 
1577  auto const& simTrack = std::visit([](auto&& obj) { return obj.g4Tracks()[0]; }, caloObj);
1578  auto const& caloPt = std::visit([](auto&& obj) { return obj.pt(); }, caloObj);
1579  stsCP_trackster_regressed_pt.push_back(caloPt);
1580 
1581  if (simTrack.crossedBoundary()) {
1582  stsCP_boundaryX.push_back(simTrack.getPositionAtBoundary().x());
1583  stsCP_boundaryY.push_back(simTrack.getPositionAtBoundary().y());
1584  stsCP_boundaryZ.push_back(simTrack.getPositionAtBoundary().z());
1585  stsCP_boundaryEta.push_back(simTrack.getPositionAtBoundary().eta());
1586  stsCP_boundaryPhi.push_back(simTrack.getPositionAtBoundary().phi());
1587  stsCP_boundaryPx.push_back(simTrack.getMomentumAtBoundary().x());
1588  stsCP_boundaryPy.push_back(simTrack.getMomentumAtBoundary().y());
1589  stsCP_boundaryPz.push_back(simTrack.getMomentumAtBoundary().z());
1590  } else {
1591  stsCP_boundaryX.push_back(-999);
1592  stsCP_boundaryY.push_back(-999);
1593  stsCP_boundaryZ.push_back(-999);
1594  stsCP_boundaryEta.push_back(-999);
1595  stsCP_boundaryPhi.push_back(-999);
1596  stsCP_boundaryPx.push_back(-999);
1597  stsCP_boundaryPy.push_back(-999);
1598  stsCP_boundaryPz.push_back(-999);
1599  }
1600  auto const trackIdx = trackster_iterator->trackIdx();
1601  stsCP_trackIdx.push_back(trackIdx);
1602  if (trackIdx != -1) {
1603  const auto& track = tracks[trackIdx];
1604 
1605  int iSide = int(track.eta() > 0);
1606 
1607  const auto& fts = trajectoryStateTransform::outerFreeState((track), bFieldProd);
1608  // to the HGCal front
1609  const auto& tsos = prop.propagate(fts, firstDisk_[iSide]->surface());
1610  if (tsos.isValid()) {
1611  const auto& globalPos = tsos.globalPosition();
1612  const auto& globalMom = tsos.globalMomentum();
1613  stsCP_track_boundaryX.push_back(globalPos.x());
1614  stsCP_track_boundaryY.push_back(globalPos.y());
1615  stsCP_track_boundaryZ.push_back(globalPos.z());
1616  stsCP_track_boundaryEta.push_back(globalPos.eta());
1617  stsCP_track_boundaryPhi.push_back(globalPos.phi());
1618  stsCP_track_boundaryPx.push_back(globalMom.x());
1619  stsCP_track_boundaryPy.push_back(globalMom.y());
1620  stsCP_track_boundaryPz.push_back(globalMom.z());
1621  stsCP_trackTime.push_back(track.t0());
1622  } else {
1623  stsCP_track_boundaryX.push_back(-999);
1624  stsCP_track_boundaryY.push_back(-999);
1625  stsCP_track_boundaryZ.push_back(-999);
1626  stsCP_track_boundaryEta.push_back(-999);
1627  stsCP_track_boundaryPhi.push_back(-999);
1628  stsCP_track_boundaryPx.push_back(-999);
1629  stsCP_track_boundaryPy.push_back(-999);
1630  stsCP_track_boundaryPz.push_back(-999);
1631  }
1632  } else {
1633  stsCP_track_boundaryX.push_back(-999);
1634  stsCP_track_boundaryY.push_back(-999);
1635  stsCP_track_boundaryZ.push_back(-999);
1636  stsCP_track_boundaryEta.push_back(-999);
1637  stsCP_track_boundaryPhi.push_back(-999);
1638  stsCP_track_boundaryPx.push_back(-999);
1639  stsCP_track_boundaryPy.push_back(-999);
1640  stsCP_track_boundaryPz.push_back(-999);
1641  }
1642  std::vector<float> id_probs;
1643  for (size_t i = 0; i < 8; i++)
1644  id_probs.push_back(trackster_iterator->id_probabilities(i));
1645  stsCP_trackster_id_probabilities.push_back(id_probs);
1646 
1647  // Clusters
1648  std::vector<uint32_t> vertices_indexes;
1649  std::vector<float> vertices_x;
1650  std::vector<float> vertices_y;
1651  std::vector<float> vertices_z;
1652  std::vector<float> vertices_time;
1653  std::vector<float> vertices_timeErr;
1654  std::vector<float> vertices_energy;
1655  std::vector<float> vertices_correctedEnergy;
1656  std::vector<float> vertices_correctedEnergyUncertainty;
1657  for (auto idx : trackster_iterator->vertices()) {
1658  vertices_indexes.push_back(idx);
1659  auto associated_cluster = (*layer_clusters_h)[idx];
1660  vertices_x.push_back(associated_cluster.x());
1661  vertices_y.push_back(associated_cluster.y());
1662  vertices_z.push_back(associated_cluster.z());
1663  vertices_energy.push_back(associated_cluster.energy());
1664  vertices_correctedEnergy.push_back(associated_cluster.correctedEnergy());
1665  vertices_correctedEnergyUncertainty.push_back(associated_cluster.correctedEnergyUncertainty());
1666  vertices_time.push_back(layerClustersTimes.get(idx).first);
1667  vertices_timeErr.push_back(layerClustersTimes.get(idx).second);
1668  }
1669  stsCP_trackster_vertices_indexes.push_back(vertices_indexes);
1670  stsCP_trackster_vertices_x.push_back(vertices_x);
1671  stsCP_trackster_vertices_y.push_back(vertices_y);
1672  stsCP_trackster_vertices_z.push_back(vertices_z);
1673  stsCP_trackster_vertices_time.push_back(vertices_time);
1674  stsCP_trackster_vertices_timeErr.push_back(vertices_timeErr);
1675  stsCP_trackster_vertices_energy.push_back(vertices_energy);
1676  stsCP_trackster_vertices_correctedEnergy.push_back(vertices_correctedEnergy);
1677  stsCP_trackster_vertices_correctedEnergyUncertainty.push_back(vertices_correctedEnergyUncertainty);
1678 
1679  // Multiplicity
1680  std::vector<float> vertices_multiplicity;
1681  for (auto multiplicity : trackster_iterator->vertex_multiplicity()) {
1682  vertices_multiplicity.push_back(multiplicity);
1683  }
1684  stsCP_trackster_vertices_multiplicity.push_back(vertices_multiplicity);
1685  }
1686 
1687  simTICLCandidate_track_in_candidate.resize(simTICLCandidates.size(), -1);
1688  for (size_t i = 0; i < simTICLCandidates.size(); ++i) {
1689  auto const& cand = simTICLCandidates[i];
1690 
1691  simTICLCandidate_raw_energy.push_back(cand.rawEnergy());
1692  simTICLCandidate_regressed_energy.push_back(cand.p4().energy());
1693  simTICLCandidate_pdgId.push_back(cand.pdgId());
1694  simTICLCandidate_charge.push_back(cand.charge());
1695  simTICLCandidate_time.push_back(cand.time());
1696  std::vector<int> tmpIdxVec;
1697  for (auto const& simTS : cand.tracksters()) {
1698  auto trackster_idx = simTS.get() - (edm::Ptr<ticl::Trackster>(simTrackstersSC_h, 0)).get();
1699  tmpIdxVec.push_back(trackster_idx);
1700  }
1701  simTICLCandidate_simTracksterCPIndex.push_back(tmpIdxVec);
1702  tmpIdxVec.clear();
1703  auto const& trackPtr = cand.trackPtr();
1704  if (!trackPtr.isNull()) {
1705  auto const& track = *trackPtr;
1706  int iSide = int(track.eta() > 0);
1707  int tk_idx = trackPtr.get() - (edm::Ptr<reco::Track>(tracks_h, 0)).get();
1709 
1710  const auto& fts = trajectoryStateTransform::outerFreeState((track), bFieldProd);
1711  // to the HGCal front
1712  const auto& tsos = prop.propagate(fts, firstDisk_[iSide]->surface());
1713  if (tsos.isValid()) {
1714  const auto& globalPos = tsos.globalPosition();
1715  const auto& globalMom = tsos.globalMomentum();
1716  simTICLCandidate_boundaryX.push_back(globalPos.x());
1717  simTICLCandidate_boundaryY.push_back(globalPos.y());
1718  simTICLCandidate_boundaryZ.push_back(globalPos.z());
1719  simTICLCandidate_boundaryPx.push_back(globalMom.x());
1720  simTICLCandidate_boundaryPy.push_back(globalMom.y());
1721  simTICLCandidate_boundaryPz.push_back(globalMom.z());
1722  } else {
1723  simTICLCandidate_boundaryX.push_back(-999);
1724  simTICLCandidate_boundaryY.push_back(-999);
1725  simTICLCandidate_boundaryZ.push_back(-999);
1726  simTICLCandidate_boundaryPx.push_back(-999);
1727  simTICLCandidate_boundaryPy.push_back(-999);
1728  simTICLCandidate_boundaryPz.push_back(-999);
1729  }
1730  } else {
1731  simTICLCandidate_boundaryX.push_back(-999);
1732  simTICLCandidate_boundaryY.push_back(-999);
1733  simTICLCandidate_boundaryZ.push_back(-999);
1734  simTICLCandidate_boundaryPx.push_back(-999);
1735  simTICLCandidate_boundaryPy.push_back(-999);
1736  simTICLCandidate_boundaryPz.push_back(-999);
1737  }
1738  }
1739 
1740  int c_id = 0;
1741 
1742  for (auto cluster_iterator = clusters.begin(); cluster_iterator != clusters.end(); ++cluster_iterator) {
1743  auto lc_seed = cluster_iterator->seed();
1744  cluster_seedID.push_back(lc_seed);
1745  cluster_energy.push_back(cluster_iterator->energy());
1746  cluster_correctedEnergy.push_back(cluster_iterator->correctedEnergy());
1747  cluster_correctedEnergyUncertainty.push_back(cluster_iterator->correctedEnergyUncertainty());
1748  cluster_position_x.push_back(cluster_iterator->x());
1749  cluster_position_y.push_back(cluster_iterator->y());
1750  cluster_position_z.push_back(cluster_iterator->z());
1751  cluster_position_eta.push_back(cluster_iterator->eta());
1752  cluster_position_phi.push_back(cluster_iterator->phi());
1753  auto haf = cluster_iterator->hitsAndFractions();
1754  auto layerId = rhtools_.getLayerWithOffset(haf[0].first);
1755  cluster_layer_id.push_back(layerId);
1756  uint32_t number_of_hits = cluster_iterator->hitsAndFractions().size();
1757  cluster_number_of_hits.push_back(number_of_hits);
1758  cluster_type.push_back(rhtools_.getCellType(lc_seed));
1759  cluster_timeErr.push_back(layerClustersTimes.get(c_id).second);
1760  cluster_time.push_back(layerClustersTimes.get(c_id).first);
1761  c_id += 1;
1762  }
1763 
1764  tracksters_in_candidate.resize(ticlcandidates.size());
1765  track_in_candidate.resize(ticlcandidates.size(), -1);
1766  nCandidates = ticlcandidates.size();
1767  for (int i = 0; i < static_cast<int>(ticlcandidates.size()); ++i) {
1768  const auto& candidate = ticlcandidates[i];
1769  candidate_charge.push_back(candidate.charge());
1770  candidate_pdgId.push_back(candidate.pdgId());
1771  candidate_energy.push_back(candidate.energy());
1772  candidate_raw_energy.push_back(candidate.rawEnergy());
1773  candidate_px.push_back(candidate.px());
1774  candidate_py.push_back(candidate.py());
1775  candidate_pz.push_back(candidate.pz());
1776  candidate_time.push_back(candidate.time());
1777  candidate_time_err.push_back(candidate.timeError());
1778  std::vector<float> id_probs;
1779  for (int j = 0; j < 8; j++) {
1781  id_probs.push_back(candidate.id_probability(type));
1782  }
1783  candidate_id_probabilities.push_back(id_probs);
1784 
1785  auto trackster_ptrs = candidate.tracksters();
1786  auto track_ptr = candidate.trackPtr();
1787  for (const auto& ts_ptr : trackster_ptrs) {
1788  auto ts_idx = ts_ptr.get() - (edm::Ptr<ticl::Trackster>(tracksters_in_candidate_handle, 0)).get();
1789  tracksters_in_candidate[i].push_back(ts_idx);
1790  }
1791  if (track_ptr.isNull())
1792  continue;
1793  int tk_idx = track_ptr.get() - (edm::Ptr<reco::Track>(tracks_h, 0)).get();
1794  track_in_candidate[i] = tk_idx;
1795  }
1796 
1797  nTrackstersMerged = trackstersmerged.size();
1798  for (auto trackster_iterator = trackstersmerged.begin(); trackster_iterator != trackstersmerged.end();
1799  ++trackster_iterator) {
1800  tracksters_merged_time.push_back(trackster_iterator->time());
1801  tracksters_merged_timeError.push_back(trackster_iterator->timeError());
1802  tracksters_merged_regressed_energy.push_back(trackster_iterator->regressed_energy());
1803  tracksters_merged_raw_energy.push_back(trackster_iterator->raw_energy());
1804  tracksters_merged_raw_em_energy.push_back(trackster_iterator->raw_em_energy());
1805  tracksters_merged_raw_pt.push_back(trackster_iterator->raw_pt());
1806  tracksters_merged_raw_em_pt.push_back(trackster_iterator->raw_em_pt());
1807  tracksters_merged_barycenter_x.push_back(trackster_iterator->barycenter().x());
1808  tracksters_merged_barycenter_y.push_back(trackster_iterator->barycenter().y());
1809  tracksters_merged_barycenter_z.push_back(trackster_iterator->barycenter().z());
1810  tracksters_merged_barycenter_eta.push_back(trackster_iterator->barycenter().eta());
1811  tracksters_merged_barycenter_phi.push_back(trackster_iterator->barycenter().phi());
1812  tracksters_merged_EV1.push_back(trackster_iterator->eigenvalues()[0]);
1813  tracksters_merged_EV2.push_back(trackster_iterator->eigenvalues()[1]);
1814  tracksters_merged_EV3.push_back(trackster_iterator->eigenvalues()[2]);
1815  tracksters_merged_eVector0_x.push_back((trackster_iterator->eigenvectors()[0]).x());
1816  tracksters_merged_eVector0_y.push_back((trackster_iterator->eigenvectors()[0]).y());
1817  tracksters_merged_eVector0_z.push_back((trackster_iterator->eigenvectors()[0]).z());
1818  tracksters_merged_sigmaPCA1.push_back(trackster_iterator->sigmasPCA()[0]);
1819  tracksters_merged_sigmaPCA2.push_back(trackster_iterator->sigmasPCA()[1]);
1820  tracksters_merged_sigmaPCA3.push_back(trackster_iterator->sigmasPCA()[2]);
1821 
1822  std::vector<float> id_probs;
1823  for (size_t i = 0; i < 8; i++)
1824  id_probs.push_back(trackster_iterator->id_probabilities(i));
1825  tracksters_merged_id_probabilities.push_back(id_probs);
1826 
1827  std::vector<uint32_t> vertices_indexes;
1828  std::vector<float> vertices_x;
1829  std::vector<float> vertices_y;
1830  std::vector<float> vertices_z;
1831  std::vector<float> vertices_time;
1832  std::vector<float> vertices_timeErr;
1833  std::vector<float> vertices_energy;
1834  std::vector<float> vertices_correctedEnergy;
1835  std::vector<float> vertices_correctedEnergyUncertainty;
1836  for (auto idx : trackster_iterator->vertices()) {
1837  vertices_indexes.push_back(idx);
1838  auto associated_cluster = (*layer_clusters_h)[idx];
1839  vertices_x.push_back(associated_cluster.x());
1840  vertices_y.push_back(associated_cluster.y());
1841  vertices_z.push_back(associated_cluster.z());
1842  vertices_energy.push_back(associated_cluster.energy());
1843  vertices_correctedEnergy.push_back(associated_cluster.correctedEnergy());
1844  vertices_correctedEnergyUncertainty.push_back(associated_cluster.correctedEnergyUncertainty());
1845  vertices_time.push_back(layerClustersTimes.get(idx).first);
1846  vertices_timeErr.push_back(layerClustersTimes.get(idx).second);
1847  }
1848  tracksters_merged_vertices_indexes.push_back(vertices_indexes);
1849  tracksters_merged_vertices_x.push_back(vertices_x);
1850  tracksters_merged_vertices_y.push_back(vertices_y);
1851  tracksters_merged_vertices_z.push_back(vertices_z);
1852  tracksters_merged_vertices_time.push_back(vertices_time);
1853  tracksters_merged_vertices_timeErr.push_back(vertices_timeErr);
1854  tracksters_merged_vertices_energy.push_back(vertices_energy);
1855  tracksters_merged_vertices_correctedEnergy.push_back(vertices_correctedEnergy);
1856  tracksters_merged_vertices_correctedEnergyUncertainty.push_back(vertices_correctedEnergyUncertainty);
1857  }
1858 
1859  // Tackster reco->sim associations
1860  trackstersCLUE3D_recoToSim_SC.resize(tracksters.size());
1861  trackstersCLUE3D_recoToSim_SC_score.resize(tracksters.size());
1862  trackstersCLUE3D_recoToSim_SC_sharedE.resize(tracksters.size());
1863  for (size_t i = 0; i < tracksters.size(); ++i) {
1864  const edm::Ref<ticl::TracksterCollection> tsRef(tracksters_handle, i);
1865 
1866  // CLUE3D -> STS-SC
1867  const auto stsSC_iter = tsRecoSimSCMap.find(tsRef);
1868  if (stsSC_iter != tsRecoSimSCMap.end()) {
1869  const auto& stsSCassociated = stsSC_iter->val;
1870  for (auto& sts : stsSCassociated) {
1871  auto sts_id = (sts.first).get() - (edm::Ref<ticl::TracksterCollection>(simTrackstersSC_h, 0)).get();
1872  trackstersCLUE3D_recoToSim_SC[i].push_back(sts_id);
1873  trackstersCLUE3D_recoToSim_SC_score[i].push_back(sts.second.second);
1874  trackstersCLUE3D_recoToSim_SC_sharedE[i].push_back(sts.second.first);
1875  }
1876  }
1877  }
1878 
1879  // SimTracksters
1880  nsimTrackstersSC = simTrackstersSC.size();
1884  for (size_t i = 0; i < nsimTrackstersSC; ++i) {
1885  const edm::Ref<ticl::TracksterCollection> stsSCRef(simTrackstersSC_h, i);
1886 
1887  // STS-SC -> CLUE3D
1888  const auto ts_iter = tsSimToRecoSCMap.find(stsSCRef);
1889  if (ts_iter != tsSimToRecoSCMap.end()) {
1890  const auto& tsAssociated = ts_iter->val;
1891  for (auto& ts : tsAssociated) {
1892  auto ts_idx = (ts.first).get() - (edm::Ref<ticl::TracksterCollection>(tracksters_handle, 0)).get();
1893  trackstersCLUE3D_simToReco_SC[i].push_back(ts_idx);
1894  trackstersCLUE3D_simToReco_SC_score[i].push_back(ts.second.second);
1895  trackstersCLUE3D_simToReco_SC_sharedE[i].push_back(ts.second.first);
1896  }
1897  }
1898  }
1899 
1900  // Tackster reco->sim associations
1901  trackstersCLUE3D_recoToSim_CP.resize(tracksters.size());
1902  trackstersCLUE3D_recoToSim_CP_score.resize(tracksters.size());
1903  trackstersCLUE3D_recoToSim_CP_sharedE.resize(tracksters.size());
1904  for (size_t i = 0; i < tracksters.size(); ++i) {
1905  const edm::Ref<ticl::TracksterCollection> tsRef(tracksters_handle, i);
1906 
1907  // CLUE3D -> STS-CP
1908  const auto stsCP_iter = tsRecoSimCPMap.find(tsRef);
1909  if (stsCP_iter != tsRecoSimCPMap.end()) {
1910  const auto& stsCPassociated = stsCP_iter->val;
1911  for (auto& sts : stsCPassociated) {
1912  auto sts_id = (sts.first).get() - (edm::Ref<ticl::TracksterCollection>(simTrackstersCP_h, 0)).get();
1913  trackstersCLUE3D_recoToSim_CP[i].push_back(sts_id);
1914  trackstersCLUE3D_recoToSim_CP_score[i].push_back(sts.second.second);
1915  trackstersCLUE3D_recoToSim_CP_sharedE[i].push_back(sts.second.first);
1916  }
1917  }
1918  }
1919 
1920  // SimTracksters
1921  nsimTrackstersCP = simTrackstersCP.size();
1925  for (size_t i = 0; i < nsimTrackstersCP; ++i) {
1926  const edm::Ref<ticl::TracksterCollection> stsCPRef(simTrackstersCP_h, i);
1927 
1928  // STS-CP -> CLUE3D
1929  const auto ts_iter = tsSimToRecoCPMap.find(stsCPRef);
1930  if (ts_iter != tsSimToRecoCPMap.end()) {
1931  const auto& tsAssociated = ts_iter->val;
1932  for (auto& ts : tsAssociated) {
1933  auto ts_idx = (ts.first).get() - (edm::Ref<ticl::TracksterCollection>(tracksters_handle, 0)).get();
1934  trackstersCLUE3D_simToReco_CP[i].push_back(ts_idx);
1935  trackstersCLUE3D_simToReco_CP_score[i].push_back(ts.second.second);
1936  trackstersCLUE3D_simToReco_CP_sharedE[i].push_back(ts.second.first);
1937  }
1938  }
1939  }
1940 
1941  // Tackster reco->sim associations
1942  MergeTracksters_recoToSim_SC.resize(trackstersmerged.size());
1943  MergeTracksters_recoToSim_SC_score.resize(trackstersmerged.size());
1944  MergeTracksters_recoToSim_SC_sharedE.resize(trackstersmerged.size());
1945  for (size_t i = 0; i < trackstersmerged.size(); ++i) {
1946  const edm::Ref<ticl::TracksterCollection> tsRef(tracksters_merged_h, i);
1947 
1948  // CLUE3D -> STS-SC
1949  const auto stsSC_iter = MergetsRecoSimSCMap.find(tsRef);
1950  if (stsSC_iter != MergetsRecoSimSCMap.end()) {
1951  const auto& stsSCassociated = stsSC_iter->val;
1952  for (auto& sts : stsSCassociated) {
1953  auto sts_id = (sts.first).get() - (edm::Ref<ticl::TracksterCollection>(simTrackstersSC_h, 0)).get();
1954  MergeTracksters_recoToSim_SC[i].push_back(sts_id);
1955  MergeTracksters_recoToSim_SC_score[i].push_back(sts.second.second);
1956  MergeTracksters_recoToSim_SC_sharedE[i].push_back(sts.second.first);
1957  }
1958  }
1959  }
1960 
1961  // SimTracksters
1962  nsimTrackstersSC = simTrackstersSC.size();
1966  for (size_t i = 0; i < nsimTrackstersSC; ++i) {
1967  const edm::Ref<ticl::TracksterCollection> stsSCRef(simTrackstersSC_h, i);
1968 
1969  // STS-SC -> CLUE3D
1970  const auto ts_iter = MergetsSimToRecoSCMap.find(stsSCRef);
1971  if (ts_iter != MergetsSimToRecoSCMap.end()) {
1972  const auto& tsAssociated = ts_iter->val;
1973  for (auto& ts : tsAssociated) {
1974  auto ts_idx = (ts.first).get() - (edm::Ref<ticl::TracksterCollection>(tracksters_merged_h, 0)).get();
1975  MergeTracksters_simToReco_SC[i].push_back(ts_idx);
1976  MergeTracksters_simToReco_SC_score[i].push_back(ts.second.second);
1977  MergeTracksters_simToReco_SC_sharedE[i].push_back(ts.second.first);
1978  }
1979  }
1980  }
1981 
1982  // Tackster reco->sim associations
1983  MergeTracksters_recoToSim_CP.resize(trackstersmerged.size());
1984  MergeTracksters_recoToSim_CP_score.resize(trackstersmerged.size());
1985  MergeTracksters_recoToSim_CP_sharedE.resize(trackstersmerged.size());
1986  for (size_t i = 0; i < trackstersmerged.size(); ++i) {
1987  const edm::Ref<ticl::TracksterCollection> tsRef(tracksters_merged_h, i);
1988 
1989  // CLUE3D -> STS-CP
1990  const auto stsCP_iter = MergetsRecoSimCPMap.find(tsRef);
1991  if (stsCP_iter != MergetsRecoSimCPMap.end()) {
1992  const auto& stsCPassociated = stsCP_iter->val;
1993  for (auto& sts : stsCPassociated) {
1994  auto sts_id = (sts.first).get() - (edm::Ref<ticl::TracksterCollection>(simTrackstersCP_h, 0)).get();
1995  MergeTracksters_recoToSim_CP[i].push_back(sts_id);
1996  MergeTracksters_recoToSim_CP_score[i].push_back(sts.second.second);
1997  MergeTracksters_recoToSim_CP_sharedE[i].push_back(sts.second.first);
1998  }
1999  }
2000  }
2001 
2002  // Tackster reco->sim associations
2003  MergeTracksters_recoToSim_PU.resize(trackstersmerged.size());
2004  MergeTracksters_recoToSim_PU_score.resize(trackstersmerged.size());
2005  MergeTracksters_recoToSim_PU_sharedE.resize(trackstersmerged.size());
2006  for (size_t i = 0; i < trackstersmerged.size(); ++i) {
2007  const edm::Ref<ticl::TracksterCollection> tsRef(tracksters_merged_h, i);
2008 
2009  // CLUE3D -> STS-PU
2010  const auto stsPU_iter = MergetsRecoSimPUMap.find(tsRef);
2011  if (stsPU_iter != MergetsRecoSimPUMap.end()) {
2012  const auto& stsPUassociated = stsPU_iter->val;
2013  for (auto& sts : stsPUassociated) {
2014  auto sts_id = (sts.first).get() - (edm::Ref<ticl::TracksterCollection>(simTrackstersPU_h, 0)).get();
2015  MergeTracksters_recoToSim_PU[i].push_back(sts_id);
2016  MergeTracksters_recoToSim_PU_score[i].push_back(sts.second.second);
2017  MergeTracksters_recoToSim_PU_sharedE[i].push_back(sts.second.first);
2018  }
2019  }
2020  }
2021 
2022  // SimTracksters
2023  nsimTrackstersCP = simTrackstersCP.size();
2027  for (size_t i = 0; i < nsimTrackstersCP; ++i) {
2028  const edm::Ref<ticl::TracksterCollection> stsCPRef(simTrackstersCP_h, i);
2029 
2030  // STS-CP -> TrackstersMerge
2031  const auto ts_iter = MergetsSimToRecoCPMap.find(stsCPRef);
2032  if (ts_iter != MergetsSimToRecoCPMap.end()) {
2033  const auto& tsAssociated = ts_iter->val;
2034  for (auto& ts : tsAssociated) {
2035  auto ts_idx = (ts.first).get() - (edm::Ref<ticl::TracksterCollection>(tracksters_merged_h, 0)).get();
2036  MergeTracksters_simToReco_CP[i].push_back(ts_idx);
2037  MergeTracksters_simToReco_CP_score[i].push_back(ts.second.second);
2038  MergeTracksters_simToReco_CP_sharedE[i].push_back(ts.second.first);
2039  }
2040  }
2041  }
2042 
2043  // SimTracksters
2044  auto nsimTrackstersPU = simTrackstersPU.size();
2045  MergeTracksters_simToReco_PU.resize(nsimTrackstersPU);
2046  MergeTracksters_simToReco_PU_score.resize(nsimTrackstersPU);
2047  MergeTracksters_simToReco_PU_sharedE.resize(nsimTrackstersPU);
2048  for (size_t i = 0; i < nsimTrackstersPU; ++i) {
2049  const edm::Ref<ticl::TracksterCollection> stsPURef(simTrackstersPU_h, i);
2050 
2051  // STS-PU -> Tracksters Merge
2052  const auto ts_iter = MergetsSimToRecoPUMap.find(stsPURef);
2053  if (ts_iter != MergetsSimToRecoPUMap.end()) {
2054  const auto& tsAssociated = ts_iter->val;
2055  for (auto& ts : tsAssociated) {
2056  auto ts_idx = (ts.first).get() - (edm::Ref<ticl::TracksterCollection>(tracksters_merged_h, 0)).get();
2057  MergeTracksters_simToReco_PU[i].push_back(ts_idx);
2058  MergeTracksters_simToReco_PU_score[i].push_back(ts.second.second);
2059  MergeTracksters_simToReco_PU_sharedE[i].push_back(ts.second.first);
2060  }
2061  }
2062  }
2063 
2064  //Tracks
2065  for (size_t i = 0; i < tracks.size(); i++) {
2066  const auto& track = tracks[i];
2067  reco::TrackRef trackref = reco::TrackRef(tracks_h, i);
2068  int iSide = int(track.eta() > 0);
2069  const auto& fts = trajectoryStateTransform::outerFreeState((track), bFieldProd);
2070  // to the HGCal front
2071  const auto& tsos = prop.propagate(fts, firstDisk_[iSide]->surface());
2072  if (tsos.isValid()) {
2073  const auto& globalPos = tsos.globalPosition();
2074  const auto& globalMom = tsos.globalMomentum();
2075  track_id.push_back(i);
2076  track_hgcal_x.push_back(globalPos.x());
2077  track_hgcal_y.push_back(globalPos.y());
2078  track_hgcal_z.push_back(globalPos.z());
2079  track_hgcal_eta.push_back(globalPos.eta());
2080  track_hgcal_phi.push_back(globalPos.phi());
2081  track_hgcal_px.push_back(globalMom.x());
2082  track_hgcal_py.push_back(globalMom.y());
2083  track_hgcal_pz.push_back(globalMom.z());
2084  track_hgcal_pt.push_back(globalMom.perp());
2085  track_pt.push_back(track.pt());
2086  track_quality.push_back(track.quality(reco::TrackBase::highPurity));
2087  track_missing_outer_hits.push_back(track.missingOuterHits());
2088  track_missing_inner_hits.push_back(track.missingInnerHits());
2089  track_charge.push_back(track.charge());
2090  track_time.push_back(trackTime[trackref]);
2091  track_time_quality.push_back(trackTimeQual[trackref]);
2092  track_time_err.push_back(trackTimeErr[trackref]);
2093  track_beta.push_back(trackBeta[trackref]);
2094  track_time_mtd.push_back(trackTimeMtd[trackref]);
2095  track_time_mtd_err.push_back(trackTimeMtdErr[trackref]);
2096  track_pos_mtd.push_back(trackPosMtd[trackref]);
2097  track_nhits.push_back(tracks[i].recHitsSize());
2098  int muId = PFMuonAlgo::muAssocToTrack(trackref, *muons_h);
2099  if (muId != -1) {
2100  const reco::MuonRef muonref = reco::MuonRef(muons_h, muId);
2101  track_isMuon.push_back(PFMuonAlgo::isMuon(muonref));
2102  track_isTrackerMuon.push_back(muons[muId].isTrackerMuon());
2103  } else {
2104  track_isMuon.push_back(-1);
2105  track_isTrackerMuon.push_back(-1);
2106  }
2107  }
2108  }
2109 
2111  trackster_tree_->Fill();
2112  if (saveLCs_)
2113  cluster_tree_->Fill();
2114  if (saveTICLCandidate_)
2115  candidate_tree_->Fill();
2117  tracksters_merged_tree_->Fill();
2118  if (saveAssociations_)
2119  associations_tree_->Fill();
2121  simtrackstersSC_tree_->Fill();
2123  simtrackstersCP_tree_->Fill();
2124  if (saveTracks_)
2125  tracks_tree_->Fill();
2127  simTICLCandidate_tree->Fill();
2128 }
2129 
2131 
2134  desc.add<edm::InputTag>("trackstersclue3d", edm::InputTag("ticlTrackstersCLUE3DHigh"));
2135  desc.add<edm::InputTag>("trackstersInCand", edm::InputTag("ticlTrackstersCLUE3DHigh"));
2136  desc.add<edm::InputTag>("layerClusters", edm::InputTag("hgcalMergeLayerClusters"));
2137  desc.add<edm::InputTag>("layer_clustersTime", edm::InputTag("hgcalMergeLayerClusters", "timeLayerCluster"));
2138  desc.add<edm::InputTag>("ticlcandidates", edm::InputTag("ticlTrackstersMerge"));
2139  desc.add<edm::InputTag>("tracks", edm::InputTag("generalTracks"));
2140  desc.add<edm::InputTag>("tracksTime", edm::InputTag("tofPID:t0"));
2141  desc.add<edm::InputTag>("tracksTimeQual", edm::InputTag("mtdTrackQualityMVA:mtdQualMVA"));
2142  desc.add<edm::InputTag>("tracksTimeErr", edm::InputTag("tofPID:sigmat0"));
2143  desc.add<edm::InputTag>("tracksBeta", edm::InputTag("trackExtenderWithMTD:generalTrackBeta"));
2144  desc.add<edm::InputTag>("tracksTimeMtd", edm::InputTag("trackExtenderWithMTD:generalTracktmtd"));
2145  desc.add<edm::InputTag>("tracksTimeMtdErr", edm::InputTag("trackExtenderWithMTD:generalTracksigmatmtd"));
2146  desc.add<edm::InputTag>("tracksPosMtd", edm::InputTag("trackExtenderWithMTD:generalTrackmtdpos"));
2147  desc.add<edm::InputTag>("trackstersmerged", edm::InputTag("ticlTrackstersMerge"));
2148  desc.add<edm::InputTag>("muons", edm::InputTag("muons1stStep"));
2149  desc.add<edm::InputTag>("simtrackstersSC", edm::InputTag("ticlSimTracksters"));
2150  desc.add<edm::InputTag>("simtrackstersCP", edm::InputTag("ticlSimTracksters", "fromCPs"));
2151  desc.add<edm::InputTag>("simtrackstersPU", edm::InputTag("ticlSimTracksters", "PU"));
2152  desc.add<edm::InputTag>("simTICLCandidates", edm::InputTag("ticlSimTracksters"));
2153  desc.add<edm::InputTag>("recoToSimAssociatorSC",
2154  edm::InputTag("tracksterSimTracksterAssociationPRbyCLUE3D", "recoToSim"));
2155  desc.add<edm::InputTag>("simToRecoAssociatorSC",
2156  edm::InputTag("tracksterSimTracksterAssociationPRbyCLUE3D", "simToReco"));
2157  desc.add<edm::InputTag>("recoToSimAssociatorCP",
2158  edm::InputTag("tracksterSimTracksterAssociationLinkingbyCLUE3D", "recoToSim"));
2159  desc.add<edm::InputTag>("simToRecoAssociatorCP",
2160  edm::InputTag("tracksterSimTracksterAssociationLinkingbyCLUE3D", "simToReco"));
2161  desc.add<edm::InputTag>("MergerecoToSimAssociatorSC",
2162  edm::InputTag("tracksterSimTracksterAssociationPR", "recoToSim"));
2163  desc.add<edm::InputTag>("MergesimToRecoAssociatorSC",
2164  edm::InputTag("tracksterSimTracksterAssociationPR", "simToReco"));
2165  desc.add<edm::InputTag>("MergerecoToSimAssociatorCP",
2166  edm::InputTag("tracksterSimTracksterAssociationLinking", "recoToSim"));
2167  desc.add<edm::InputTag>("MergesimToRecoAssociatorCP",
2168  edm::InputTag("tracksterSimTracksterAssociationLinking", "simToReco"));
2169  desc.add<edm::InputTag>("MergerecoToSimAssociatorPU",
2170  edm::InputTag("tracksterSimTracksterAssociationLinkingPU", "recoToSim"));
2171  desc.add<edm::InputTag>("MergesimToRecoAssociatorPU",
2172  edm::InputTag("tracksterSimTracksterAssociationLinkingPU", "simToReco"));
2173  desc.add<edm::InputTag>("simclusters", edm::InputTag("mix", "MergedCaloTruth"));
2174  desc.add<edm::InputTag>("caloparticles", edm::InputTag("mix", "MergedCaloTruth"));
2175  desc.add<std::string>("detector", "HGCAL");
2176  desc.add<std::string>("propagator", "PropagatorWithMaterial");
2177 
2178  desc.add<bool>("saveLCs", true);
2179  desc.add<bool>("saveCLUE3DTracksters", true);
2180  desc.add<bool>("saveTrackstersMerged", true);
2181  desc.add<bool>("saveSimTrackstersSC", true);
2182  desc.add<bool>("saveSimTrackstersCP", true);
2183  desc.add<bool>("saveTICLCandidate", true);
2184  desc.add<bool>("saveSimTICLCandidate", true);
2185  desc.add<bool>("saveTracks", true);
2186  desc.add<bool>("saveAssociations", true);
2187  descriptions.add("ticlDumper", desc);
2188 }
2189 
std::vector< float > stsSC_trackster_time
Definition: TICLDumper.cc:195
size_t nsimTrackstersCP
Definition: TICLDumper.cc:160
TTree * simtrackstersSC_tree_
Definition: TICLDumper.cc:449
std::vector< std::vector< float > > MergeTracksters_recoToSim_PU_score
Definition: TICLDumper.cc:397
std::vector< float > stsCP_trackster_EV3
Definition: TICLDumper.cc:264
double waferZ(int layer, bool reco) const
std::vector< float > stsSC_trackster_sigmaPCA1
Definition: TICLDumper.cc:215
const edm::EDGetTokenT< std::vector< ticl::Trackster > > simTracksters_PU_token_
Definition: TICLDumper.cc:111
std::vector< float > stsSC_trackster_barycenter_x
Definition: TICLDumper.cc:204
std::vector< float > stsCP_trackster_eVector0_x
Definition: TICLDumper.cc:265
std::vector< float > cluster_timeErr
Definition: TICLDumper.cc:415
static DiskPointer build(Args &&... args)
Definition: BoundDisk.h:38
const edm::EDGetTokenT< std::vector< double > > hgcaltracks_eta_token_
Definition: TICLDumper.cc:99
std::vector< std::vector< float > > MergeTracksters_simToReco_PU_score
Definition: TICLDumper.cc:400
const edm::ESGetToken< CaloGeometry, CaloGeometryRecord > geometry_token_
Definition: TICLDumper.cc:126
static int muAssocToTrack(const reco::TrackRef &trackref, const reco::MuonCollection &muons)
Definition: PFMuonAlgo.cc:479
std::vector< float > cluster_position_x
Definition: TICLDumper.cc:407
edm::AssociationMap< edm::OneToManyWithQualityGeneric< ticl::TracksterCollection, ticl::TracksterCollection, std::pair< float, float > > > RecoToSimCollectionSimTracksters
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
std::vector< float > stsCP_track_boundaryPx
Definition: TICLDumper.cc:287
std::vector< int > simTICLCandidate_charge
Definition: TICLDumper.cc:314
std::vector< std::vector< float > > stsSC_trackster_vertices_correctedEnergy
Definition: TICLDumper.cc:245
std::vector< std::vector< float > > tracksters_merged_vertices_energy
Definition: TICLDumper.cc:361
std::vector< float > stsCP_trackster_eVector0_z
Definition: TICLDumper.cc:267
TTree * cluster_tree_
Definition: TICLDumper.cc:445
std::vector< float > stsCP_trackster_barycenter_eta
Definition: TICLDumper.cc:260
std::vector< float > stsSC_track_boundaryPz
Definition: TICLDumper.cc:236
std::vector< float > trackster_EV1
Definition: TICLDumper.cc:172
std::vector< double > candidate_py
Definition: TICLDumper.cc:324
std::vector< float > stsSC_boundaryPx
Definition: TICLDumper.cc:226
std::vector< float > tracksters_merged_eVector0_x
Definition: TICLDumper.cc:349
ProductID id() const
Definition: HandleBase.cc:29
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
std::vector< float > stsCP_trackster_timeBoundary
Definition: TICLDumper.cc:249
const edm::EDGetTokenT< std::vector< double > > hgcaltracks_py_token_
Definition: TICLDumper.cc:102
std::vector< float > trackster_sigmaPCA1
Definition: TICLDumper.cc:178
bool saveTrackstersMerged_
Definition: TICLDumper.cc:140
std::vector< float > trackster_barycenter_y
Definition: TICLDumper.cc:170
unsigned int stsCP_ntracksters_
Definition: TICLDumper.cc:158
bool saveAssociations_
Definition: TICLDumper.cc:146
std::vector< float > simTICLCandidate_raw_energy
Definition: TICLDumper.cc:302
void endRun(edm::Run const &iEvent, edm::EventSetup const &) override
Definition: TICLDumper.cc:79
std::vector< float > stsSC_trackster_barycenter_y
Definition: TICLDumper.cc:205
TTree * simtrackstersCP_tree_
Definition: TICLDumper.cc:450
const edm::EDGetTokenT< edm::ValueMap< float > > tracks_beta_token_
Definition: TICLDumper.cc:92
std::vector< std::vector< float > > trackster_vertices_correctedEnergyUncertainty
Definition: TICLDumper.cc:192
std::vector< std::vector< float > > trackster_vertices_multiplicity
Definition: TICLDumper.cc:193
const edm::EDGetTokenT< edm::ValueMap< GlobalPoint > > tracks_pos_mtd_token_
Definition: TICLDumper.cc:95
static bool isMuon(const reco::PFBlockElement &elt)
Definition: PFMuonAlgo.cc:48
std::vector< float > stsSC_boundaryX
Definition: TICLDumper.cc:221
std::vector< float > stsSC_track_boundaryPhi
Definition: TICLDumper.cc:233
std::vector< float > stsSC_trackster_timeError
Definition: TICLDumper.cc:197
std::vector< std::vector< float > > trackster_vertices_y
Definition: TICLDumper.cc:186
bool saveSimTrackstersSC_
Definition: TICLDumper.cc:141
std::vector< float > stsCP_boundaryY
Definition: TICLDumper.cc:275
std::vector< std::vector< uint32_t > > tracksters_in_candidate
Definition: TICLDumper.cc:329
std::vector< float > cluster_position_y
Definition: TICLDumper.cc:408
std::vector< int > stsSC_pdgID
Definition: TICLDumper.cc:218
std::vector< float > stsCP_boundaryX
Definition: TICLDumper.cc:274
std::vector< std::vector< float > > stsCP_trackster_vertices_time
Definition: TICLDumper.cc:295
void endJob() override
Definition: TICLDumper.cc:2130
std::vector< float > cluster_correctedEnergy
Definition: TICLDumper.cc:405
std::vector< float > stsSC_trackster_regressed_pt
Definition: TICLDumper.cc:199
std::vector< float > track_beta
Definition: TICLDumper.cc:436
std::vector< float > trackster_raw_pt
Definition: TICLDumper.cc:167
std::vector< float > trackster_eVector0_z
Definition: TICLDumper.cc:177
std::vector< uint32_t > cluster_seedID
Definition: TICLDumper.cc:403
void buildLayers()
Definition: TICLDumper.cc:1144
std::vector< std::vector< uint32_t > > MergeTracksters_recoToSim_SC
Definition: TICLDumper.cc:382
std::vector< double > Vec
Definition: TICLDumper.cc:66
std::vector< std::vector< float > > MergeTracksters_simToReco_SC_sharedE
Definition: TICLDumper.cc:387
int getCellType(const DetId &id) const
Definition: RecHitTools.cc:444
std::vector< float > cluster_position_z
Definition: TICLDumper.cc:409
std::vector< std::vector< float > > stsSC_trackster_vertices_z
Definition: TICLDumper.cc:241
std::vector< int > track_quality
Definition: TICLDumper.cc:429
std::vector< float > stsSC_trackster_eVector0_z
Definition: TICLDumper.cc:214
std::vector< float > track_hgcal_y
Definition: TICLDumper.cc:420
std::vector< std::vector< float > > MergeTracksters_recoToSim_SC_score
Definition: TICLDumper.cc:383
std::vector< std::vector< uint32_t > > stsSC_trackster_vertices_indexes
Definition: TICLDumper.cc:238
size_t nsimTrackstersSC
Definition: TICLDumper.cc:159
std::vector< unsigned int > track_id
Definition: TICLDumper.cc:418
const edm::EDGetTokenT< std::vector< TICLCandidate > > ticl_candidates_token_
Definition: TICLDumper.cc:86
std::vector< float > stsCP_trackster_barycenter_phi
Definition: TICLDumper.cc:261
std::vector< float > tracksters_merged_sigmaPCA3
Definition: TICLDumper.cc:354
std::vector< float > simTICLCandidate_boundaryY
Definition: TICLDumper.cc:306
std::vector< std::vector< float > > trackster_vertices_energy
Definition: TICLDumper.cc:190
std::vector< float > trackster_raw_em_energy
Definition: TICLDumper.cc:166
const edm::EDGetTokenT< std::vector< ticl::Trackster > > simTracksters_SC_token_
Definition: TICLDumper.cc:109
std::vector< std::vector< uint32_t > > trackstersCLUE3D_simToReco_CP
Definition: TICLDumper.cc:378
const edm::EDGetTokenT< edm::ValueMap< float > > tracks_time_quality_token_
Definition: TICLDumper.cc:90
std::vector< std::vector< float > > trackster_vertices_time
Definition: TICLDumper.cc:188
std::vector< float > stsSC_boundaryPz
Definition: TICLDumper.cc:228
std::vector< std::vector< float > > trackster_vertices_z
Definition: TICLDumper.cc:187
std::vector< float > stsCP_trackster_EV1
Definition: TICLDumper.cc:262
void clearVariables()
Definition: TICLDumper.cc:455
const HGCalDDDConstants * hgcons_
Definition: TICLDumper.cc:133
std::vector< float > cluster_energy
Definition: TICLDumper.cc:404
const edm::ESGetToken< Propagator, TrackingComponentsRecord > propagator_token_
Definition: TICLDumper.cc:130
std::vector< std::vector< float > > MergeTracksters_recoToSim_CP_sharedE
Definition: TICLDumper.cc:391
std::vector< int > candidate_charge
Definition: TICLDumper.cc:319
std::vector< float > track_hgcal_pz
Definition: TICLDumper.cc:424
const edm::EDGetTokenT< ticl::RecoToSimCollectionSimTracksters > tsRecoToSimSC_token_
Definition: TICLDumper.cc:113
std::vector< float > stsSC_trackster_timeBoundary
Definition: TICLDumper.cc:196
std::vector< float > trackster_barycenter_x
Definition: TICLDumper.cc:169
std::unique_ptr< GeomDet > interfaceDisk_[2]
Definition: TICLDumper.cc:135
muons
the two sets of parameters below are mutually exclusive, depending if RECO or ALCARECO is used the us...
Definition: DiMuonV_cfg.py:214
std::vector< std::vector< float > > stsSC_trackster_vertices_time
Definition: TICLDumper.cc:242
TrajectoryStateOnSurface propagate(STA const &state, SUR const &surface) const
Definition: Propagator.h:50
void beginJob() override
Definition: TICLDumper.cc:829
unsigned int stsSC_ntracksters_
Definition: TICLDumper.cc:157
std::vector< std::vector< float > > trackster_id_probabilities
Definition: TICLDumper.cc:183
std::vector< float > stsCP_trackTime
Definition: TICLDumper.cc:273
std::vector< float > stsCP_boundaryEta
Definition: TICLDumper.cc:277
std::vector< std::vector< float > > tracksters_merged_vertices_z
Definition: TICLDumper.cc:358
std::vector< std::vector< float > > stsCP_trackster_vertices_timeErr
Definition: TICLDumper.cc:296
const edm::EDGetTokenT< ticl::SimToRecoCollectionSimTracksters > tsSimToRecoCP_token_
Definition: TICLDumper.cc:116
std::vector< std::vector< uint32_t > > MergeTracksters_recoToSim_PU
Definition: TICLDumper.cc:396
std::vector< float > candidate_time
Definition: TICLDumper.cc:326
std::vector< float > stsSC_trackster_raw_pt
Definition: TICLDumper.cc:202
const edm::EDGetTokenT< ticl::RecoToSimCollectionSimTracksters > tsRecoToSimCP_token_
Definition: TICLDumper.cc:115
bool saveSimTrackstersCP_
Definition: TICLDumper.cc:142
std::vector< std::vector< uint32_t > > trackstersCLUE3D_recoToSim_SC
Definition: TICLDumper.cc:368
std::vector< float > tracksters_merged_time
Definition: TICLDumper.cc:334
std::vector< float > stsCP_boundaryPz
Definition: TICLDumper.cc:281
std::vector< float > tracksters_merged_sigmaPCA2
Definition: TICLDumper.cc:353
std::vector< float > track_hgcal_px
Definition: TICLDumper.cc:422
TTree * candidate_tree_
Definition: TICLDumper.cc:446
std::vector< std::vector< float > > MergeTracksters_recoToSim_SC_sharedE
Definition: TICLDumper.cc:384
std::vector< unsigned int > cluster_layer_id
Definition: TICLDumper.cc:412
std::vector< float > candidate_energy
Definition: TICLDumper.cc:321
std::vector< float > tracksters_merged_barycenter_phi
Definition: TICLDumper.cc:345
std::vector< int > track_charge
Definition: TICLDumper.cc:432
std::vector< float > stsSC_boundaryPhi
Definition: TICLDumper.cc:225
std::vector< float > tracksters_merged_eVector0_z
Definition: TICLDumper.cc:351
std::vector< float > stsCP_track_boundaryPy
Definition: TICLDumper.cc:288
std::vector< float > tracksters_merged_raw_pt
Definition: TICLDumper.cc:339
std::vector< float > tracksters_merged_raw_em_energy
Definition: TICLDumper.cc:338
std::vector< int > simTICLCandidate_track_in_candidate
Definition: TICLDumper.cc:315
std::vector< float > track_hgcal_eta
Definition: TICLDumper.cc:425
std::vector< float > tracksters_merged_raw_em_pt
Definition: TICLDumper.cc:340
std::vector< float > stsSC_track_boundaryPx
Definition: TICLDumper.cc:234
std::vector< int > cluster_type
Definition: TICLDumper.cc:413
std::vector< std::vector< float > > tracksters_merged_vertices_multiplicity
Definition: TICLDumper.cc:364
std::vector< std::vector< uint32_t > > stsCP_trackster_vertices_indexes
Definition: TICLDumper.cc:291
const edm::EDGetTokenT< std::vector< CaloParticle > > caloparticles_token_
Definition: TICLDumper.cc:124
std::vector< float > trackster_barycenter_z
Definition: TICLDumper.cc:171
std::vector< float > stsCP_trackster_sigmaPCA1
Definition: TICLDumper.cc:268
const edm::EDGetTokenT< edm::ValueMap< float > > tracks_time_token_
Definition: TICLDumper.cc:89
size_t nCandidates
Definition: TICLDumper.cc:318
int iEvent
Definition: GenABIO.cc:224
T const * product() const
Definition: ESHandle.h:86
std::vector< std::vector< float > > stsCP_trackster_vertices_y
Definition: TICLDumper.cc:293
std::vector< float > simTICLCandidate_boundaryZ
Definition: TICLDumper.cc:307
std::vector< float > stsCP_trackster_barycenter_z
Definition: TICLDumper.cc:259
const edm::EDGetTokenT< std::vector< ticl::Trackster > > tracksters_in_candidate_token_
Definition: TICLDumper.cc:84
const edm::EDGetTokenT< std::vector< ticl::Trackster > > tracksters_token_
Definition: TICLDumper.cc:83
GlobalPoint globalPosition() const
std::vector< std::vector< uint32_t > > MergeTracksters_simToReco_PU
Definition: TICLDumper.cc:399
std::vector< float > track_hgcal_z
Definition: TICLDumper.cc:421
std::vector< std::vector< float > > trackstersCLUE3D_simToReco_SC_sharedE
Definition: TICLDumper.cc:373
std::vector< std::vector< float > > stsCP_trackster_vertices_correctedEnergyUncertainty
Definition: TICLDumper.cc:299
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: TICLDumper.cc:1181
std::vector< float > stsCP_boundaryZ
Definition: TICLDumper.cc:276
Definition: Muon.py:1
std::vector< float > stsCP_track_boundaryPhi
Definition: TICLDumper.cc:286
std::vector< std::vector< float > > stsCP_trackster_vertices_z
Definition: TICLDumper.cc:294
std::vector< int > track_missing_outer_hits
Definition: TICLDumper.cc:430
std::vector< float > stsSC_trackster_raw_em_energy
Definition: TICLDumper.cc:201
std::vector< float > tracksters_merged_barycenter_eta
Definition: TICLDumper.cc:344
std::vector< float > stsSC_track_boundaryX
Definition: TICLDumper.cc:229
std::vector< std::vector< float > > MergeTracksters_recoToSim_CP_score
Definition: TICLDumper.cc:390
std::vector< float > tracksters_merged_barycenter_y
Definition: TICLDumper.cc:342
std::vector< float > track_hgcal_phi
Definition: TICLDumper.cc:426
std::vector< std::vector< float > > tracksters_merged_vertices_timeErr
Definition: TICLDumper.cc:360
std::vector< float > cluster_time
Definition: TICLDumper.cc:414
std::vector< float > simTICLCandidate_boundaryPy
Definition: TICLDumper.cc:309
std::vector< double > track_time
Definition: TICLDumper.cc:433
TTree * tracksters_merged_tree_
Definition: TICLDumper.cc:447
const edm::EDGetTokenT< ticl::SimToRecoCollectionSimTracksters > MergeSimToRecoCP_token_
Definition: TICLDumper.cc:120
const edm::EDGetTokenT< ticl::RecoToSimCollectionSimTracksters > MergeRecoToSimPU_token_
Definition: TICLDumper.cc:121
std::vector< double > candidate_pz
Definition: TICLDumper.cc:325
std::vector< std::vector< float > > stsSC_trackster_vertices_y
Definition: TICLDumper.cc:240
std::vector< std::vector< float > > MergeTracksters_simToReco_CP_sharedE
Definition: TICLDumper.cc:394
TTree * associations_tree_
Definition: TICLDumper.cc:448
std::vector< float > trackster_barycenter_phi
Definition: TICLDumper.cc:182
std::vector< float > stsCP_trackster_timeError
Definition: TICLDumper.cc:250
std::vector< std::vector< float > > tracksters_merged_vertices_correctedEnergyUncertainty
Definition: TICLDumper.cc:363
std::vector< std::vector< float > > trackstersCLUE3D_recoToSim_SC_sharedE
Definition: TICLDumper.cc:370
std::vector< float > stsSC_trackTime
Definition: TICLDumper.cc:220
std::vector< std::vector< float > > stsCP_trackster_id_probabilities
Definition: TICLDumper.cc:290
Transition
Definition: Transition.h:12
bool saveTICLCandidate_
Definition: TICLDumper.cc:143
std::vector< float > stsSC_track_boundaryZ
Definition: TICLDumper.cc:231
std::vector< float > trackster_sigmaPCA2
Definition: TICLDumper.cc:179
edm::Ref< MuonCollection > MuonRef
presistent reference to a Muon
Definition: MuonFwd.h:13
std::vector< float > tracksters_merged_sigmaPCA1
Definition: TICLDumper.cc:352
std::vector< std::vector< float > > MergeTracksters_recoToSim_PU_sharedE
Definition: TICLDumper.cc:398
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::ESHandle< Propagator > propagator_
Definition: TICLDumper.cc:137
std::vector< float > simTICLCandidate_boundaryPx
Definition: TICLDumper.cc:308
std::vector< float > stsSC_trackster_raw_energy
Definition: TICLDumper.cc:200
std::vector< float > simTICLCandidate_boundaryX
Definition: TICLDumper.cc:305
std::vector< double > candidate_px
Definition: TICLDumper.cc:323
std::vector< float > stsCP_trackster_EV2
Definition: TICLDumper.cc:263
const edm::EDGetTokenT< std::vector< double > > hgcaltracks_x_token_
Definition: TICLDumper.cc:96
std::vector< float > tracksters_merged_EV1
Definition: TICLDumper.cc:346
std::vector< float > tracksters_merged_timeError
Definition: TICLDumper.cc:335
std::vector< float > stsSC_trackster_EV3
Definition: TICLDumper.cc:211
std::vector< float > stsCP_track_boundaryZ
Definition: TICLDumper.cc:284
std::vector< std::vector< uint32_t > > trackster_vertices_indexes
Definition: TICLDumper.cc:184
size_t nTrackstersMerged
Definition: TICLDumper.cc:333
const edm::EDGetTokenT< std::vector< double > > hgcaltracks_z_token_
Definition: TICLDumper.cc:98
std::vector< float > track_hgcal_pt
Definition: TICLDumper.cc:427
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
std::vector< int > simTICLCandidate_pdgId
Definition: TICLDumper.cc:313
TICLDumper(const edm::ParameterSet &)
Definition: TICLDumper.cc:745
const edm::EDGetTokenT< std::vector< double > > hgcaltracks_px_token_
Definition: TICLDumper.cc:101
std::vector< std::vector< float > > stsSC_trackster_vertices_x
Definition: TICLDumper.cc:239
std::vector< float > stsCP_boundaryPhi
Definition: TICLDumper.cc:278
std::vector< float > trackster_raw_em_pt
Definition: TICLDumper.cc:168
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: TICLDumper.cc:2132
std::vector< float > stsSC_track_boundaryPy
Definition: TICLDumper.cc:235
std::vector< std::vector< uint32_t > > MergeTracksters_simToReco_SC
Definition: TICLDumper.cc:385
std::vector< float > stsSC_boundaryPy
Definition: TICLDumper.cc:227
std::vector< float > stsCP_trackster_time
Definition: TICLDumper.cc:248
std::vector< std::vector< float > > trackstersCLUE3D_recoToSim_CP_sharedE
Definition: TICLDumper.cc:377
std::vector< float > track_time_quality
Definition: TICLDumper.cc:434
const std::string propName_
Definition: TICLDumper.cc:128
std::pair< double, double > rangeR(double z, bool reco) const
std::vector< float > tracksters_merged_regressed_energy
Definition: TICLDumper.cc:336
edm::AssociationMap< edm::OneToManyWithQualityGeneric< ticl::TracksterCollection, ticl::TracksterCollection, std::pair< float, float > > > SimToRecoCollectionSimTracksters
std::vector< std::vector< float > > tracksters_merged_vertices_y
Definition: TICLDumper.cc:357
bool saveLCs_
Definition: TICLDumper.cc:138
std::vector< float > trackster_eVector0_x
Definition: TICLDumper.cc:175
TTree * tracks_tree_
Definition: TICLDumper.cc:451
std::vector< float > trackster_eVector0_y
Definition: TICLDumper.cc:176
std::vector< float > trackster_EV3
Definition: TICLDumper.cc:174
std::vector< float > stsSC_trackster_barycenter_eta
Definition: TICLDumper.cc:207
std::vector< float > trackster_time
Definition: TICLDumper.cc:162
const edm::EDGetTokenT< ticl::RecoToSimCollectionSimTracksters > MergeRecoToSimCP_token_
Definition: TICLDumper.cc:119
std::vector< std::vector< float > > stsCP_trackster_vertices_correctedEnergy
Definition: TICLDumper.cc:298
std::vector< float > simTICLCandidate_regressed_energy
Definition: TICLDumper.cc:303
std::vector< std::vector< float > > trackstersCLUE3D_simToReco_CP_score
Definition: TICLDumper.cc:379
std::vector< std::vector< float > > tracksters_merged_vertices_time
Definition: TICLDumper.cc:359
const edm::EDGetTokenT< std::vector< bool > > tracks_mask_token_
Definition: TICLDumper.cc:88
std::vector< std::vector< float > > stsSC_trackster_vertices_correctedEnergyUncertainty
Definition: TICLDumper.cc:246
std::vector< float > simTICLCandidate_caloParticleMass
Definition: TICLDumper.cc:311
std::vector< float > stsSC_trackster_barycenter_z
Definition: TICLDumper.cc:206
std::vector< std::vector< float > > stsSC_trackster_vertices_energy
Definition: TICLDumper.cc:244
std::vector< float > stsSC_track_boundaryY
Definition: TICLDumper.cc:230
std::vector< std::vector< float > > trackstersCLUE3D_simToReco_SC_score
Definition: TICLDumper.cc:372
const edm::EDGetTokenT< std::vector< int > > tracksterSeeds_token_
Definition: TICLDumper.cc:107
std::vector< std::vector< float > > tracksters_merged_id_probabilities
Definition: TICLDumper.cc:365
void beginRun(const edm::Run &, const edm::EventSetup &) override
Definition: TICLDumper.cc:817
edm::ESGetToken< HGCalDDDConstants, IdealGeometryRecord > hdc_token_
Definition: TICLDumper.cc:132
std::vector< float > stsCP_trackster_sigmaPCA3
Definition: TICLDumper.cc:270
bool saveTracks_
Definition: TICLDumper.cc:145
const edm::EDGetTokenT< std::vector< double > > hgcaltracks_y_token_
Definition: TICLDumper.cc:97
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:20
std::vector< float > tracksters_merged_EV2
Definition: TICLDumper.cc:347
const edm::EDGetTokenT< edm::ValueMap< float > > tracks_time_mtd_token_
Definition: TICLDumper.cc:93
std::vector< std::vector< uint32_t > > MergeTracksters_recoToSim_CP
Definition: TICLDumper.cc:389
std::vector< float > stsSC_trackster_sigmaPCA3
Definition: TICLDumper.cc:217
const edm::EDGetTokenT< std::vector< ticl::Trackster > > simTracksters_CP_token_
Definition: TICLDumper.cc:110
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< std::vector< float > > stsSC_trackster_vertices_multiplicity
Definition: TICLDumper.cc:247
GlobalPoint getPositionLayer(int layer, bool nose=false) const
Definition: RecHitTools.cc:152
std::vector< float > stsSC_boundaryZ
Definition: TICLDumper.cc:223
std::vector< std::vector< float > > stsSC_trackster_id_probabilities
Definition: TICLDumper.cc:237
std::vector< float > stsSC_trackster_EV1
Definition: TICLDumper.cc:209
std::vector< float > cluster_position_eta
Definition: TICLDumper.cc:410
const edm::EDGetTokenT< std::vector< ticl::Trackster > > tracksters_merged_token_
Definition: TICLDumper.cc:104
hgcal::RecHitTools rhtools_
Definition: TICLDumper.cc:131
std::vector< float > track_pt
Definition: TICLDumper.cc:428
unsigned int nclusters_
Definition: TICLDumper.cc:156
std::vector< std::vector< float > > candidate_id_probabilities
Definition: TICLDumper.cc:328
std::vector< int > track_isTrackerMuon
Definition: TICLDumper.cc:442
std::vector< float > trackster_barycenter_eta
Definition: TICLDumper.cc:181
std::vector< float > stsCP_trackster_regressed_energy
Definition: TICLDumper.cc:251
std::vector< float > stsCP_trackster_sigmaPCA2
Definition: TICLDumper.cc:269
std::vector< float > stsCP_track_boundaryPz
Definition: TICLDumper.cc:289
void setGeometry(CaloGeometry const &)
Definition: RecHitTools.cc:79
ticl::Vector Vector
Definition: TICLDumper.cc:65
std::vector< std::vector< float > > MergeTracksters_simToReco_PU_sharedE
Definition: TICLDumper.cc:401
std::vector< float > stsSC_trackster_regressed_energy
Definition: TICLDumper.cc:198
std::vector< std::vector< uint32_t > > MergeTracksters_simToReco_CP
Definition: TICLDumper.cc:392
std::vector< std::vector< int > > simTICLCandidate_simTracksterCPIndex
Definition: TICLDumper.cc:304
std::vector< float > tracksters_merged_EV3
Definition: TICLDumper.cc:348
std::vector< float > track_hgcal_py
Definition: TICLDumper.cc:423
bool saveSimTICLCandidate_
Definition: TICLDumper.cc:144
std::vector< int > track_in_candidate
Definition: TICLDumper.cc:330
math::XYZVectorF Vector
Definition: Common.h:42
std::vector< float > track_time_err
Definition: TICLDumper.cc:435
std::vector< std::vector< float > > trackster_vertices_correctedEnergy
Definition: TICLDumper.cc:191
std::vector< int > track_isMuon
Definition: TICLDumper.cc:441
fixed size matrix
std::vector< std::vector< uint32_t > > trackstersCLUE3D_simToReco_SC
Definition: TICLDumper.cc:371
HLT enums.
std::vector< float > stsSC_trackster_raw_em_pt
Definition: TICLDumper.cc:203
std::vector< std::vector< float > > trackstersCLUE3D_recoToSim_SC_score
Definition: TICLDumper.cc:369
std::vector< float > trackster_sigmaPCA3
Definition: TICLDumper.cc:180
std::vector< int > stsCP_trackIdx
Definition: TICLDumper.cc:272
std::vector< float > candidate_raw_energy
Definition: TICLDumper.cc:322
const edm::EDGetTokenT< edm::ValueMap< float > > tracks_time_mtd_err_token_
Definition: TICLDumper.cc:94
std::vector< float > cluster_correctedEnergyUncertainty
Definition: TICLDumper.cc:406
std::vector< float > stsCP_trackster_barycenter_x
Definition: TICLDumper.cc:257
std::vector< int > candidate_pdgId
Definition: TICLDumper.cc:320
std::unique_ptr< GeomDet > firstDisk_[2]
Definition: TICLDumper.cc:134
std::vector< float > stsCP_trackster_regressed_pt
Definition: TICLDumper.cc:252
std::vector< std::vector< float > > MergeTracksters_simToReco_CP_score
Definition: TICLDumper.cc:393
std::vector< float > stsCP_trackster_raw_em_energy
Definition: TICLDumper.cc:254
std::vector< int > stsCP_pdgID
Definition: TICLDumper.cc:271
Definition: Common.h:10
const edm::EDGetTokenT< ticl::SimToRecoCollectionSimTracksters > MergeSimToRecoPU_token_
Definition: TICLDumper.cc:122
std::vector< float > track_time_mtd_err
Definition: TICLDumper.cc:438
std::vector< float > stsCP_track_boundaryX
Definition: TICLDumper.cc:282
std::vector< std::vector< float > > stsCP_trackster_vertices_x
Definition: TICLDumper.cc:292
std::vector< std::vector< float > > MergeTracksters_simToReco_SC_score
Definition: TICLDumper.cc:386
const edm::EDGetTokenT< edm::ValueMap< std::pair< float, float > > > clustersTime_token_
Definition: TICLDumper.cc:106
const edm::EDGetTokenT< edm::ValueMap< float > > tracks_time_err_token_
Definition: TICLDumper.cc:91
std::vector< float > trackster_regressed_energy
Definition: TICLDumper.cc:164
std::vector< std::vector< uint32_t > > tracksters_merged_vertices_indexes
Definition: TICLDumper.cc:355
std::vector< std::vector< uint32_t > > trackstersCLUE3D_recoToSim_CP
Definition: TICLDumper.cc:375
std::vector< std::vector< float > > trackstersCLUE3D_simToReco_CP_sharedE
Definition: TICLDumper.cc:380
const edm::EDGetTokenT< std::vector< double > > hgcaltracks_phi_token_
Definition: TICLDumper.cc:100
std::vector< int > track_nhits
Definition: TICLDumper.cc:440
const edm::EDGetTokenT< std::vector< TICLCandidate > > simTICLCandidate_token_
Definition: TICLDumper.cc:112
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > bfield_token_
Definition: TICLDumper.cc:129
std::vector< int > track_missing_inner_hits
Definition: TICLDumper.cc:431
std::vector< std::vector< float > > trackster_vertices_x
Definition: TICLDumper.cc:185
std::vector< float > track_hgcal_x
Definition: TICLDumper.cc:419
std::vector< GlobalPoint > track_pos_mtd
Definition: TICLDumper.cc:439
const edm::EDGetTokenT< ticl::SimToRecoCollectionSimTracksters > tsSimToRecoSC_token_
Definition: TICLDumper.cc:114
std::vector< float > stsSC_boundaryY
Definition: TICLDumper.cc:222
std::vector< float > track_time_mtd
Definition: TICLDumper.cc:437
std::vector< float > stsCP_trackster_barycenter_y
Definition: TICLDumper.cc:258
std::vector< float > stsCP_boundaryPx
Definition: TICLDumper.cc:279
std::vector< float > stsSC_trackster_eVector0_x
Definition: TICLDumper.cc:212
std::vector< float > stsCP_trackster_raw_energy
Definition: TICLDumper.cc:253
std::vector< float > tracksters_merged_raw_energy
Definition: TICLDumper.cc:337
std::vector< float > stsSC_track_boundaryEta
Definition: TICLDumper.cc:232
std::vector< uint32_t > cluster_number_of_hits
Definition: TICLDumper.cc:416
std::vector< std::vector< float > > trackstersCLUE3D_recoToSim_CP_score
Definition: TICLDumper.cc:376
TTree * simTICLCandidate_tree
Definition: TICLDumper.cc:452
const edm::EDGetTokenT< ticl::SimToRecoCollectionSimTracksters > MergeSimToRecoSC_token_
Definition: TICLDumper.cc:118
std::vector< std::vector< float > > trackster_vertices_timeErr
Definition: TICLDumper.cc:189
std::vector< float > stsCP_trackster_raw_em_pt
Definition: TICLDumper.cc:256
void initialize(const HGCalDDDConstants *hgcons, const hgcal::RecHitTools rhtools, const edm::ESHandle< MagneticField > bfieldH, const edm::ESHandle< Propagator > propH)
Definition: TICLDumper.cc:1170
std::vector< float > stsCP_trackster_eVector0_y
Definition: TICLDumper.cc:266
const edm::EDGetTokenT< std::vector< reco::CaloCluster > > layer_clusters_token_
Definition: TICLDumper.cc:85
std::vector< float > stsSC_trackster_EV2
Definition: TICLDumper.cc:210
std::vector< float > candidate_time_err
Definition: TICLDumper.cc:327
std::vector< float > trackster_raw_energy
Definition: TICLDumper.cc:165
std::vector< float > stsSC_trackster_eVector0_y
Definition: TICLDumper.cc:213
std::vector< float > trackster_EV2
Definition: TICLDumper.cc:173
std::vector< float > tracksters_merged_eVector0_y
Definition: TICLDumper.cc:350
std::vector< std::vector< float > > stsCP_trackster_vertices_multiplicity
Definition: TICLDumper.cc:300
const edm::EDGetTokenT< std::vector< double > > hgcaltracks_pz_token_
Definition: TICLDumper.cc:103
std::vector< std::vector< float > > stsSC_trackster_vertices_timeErr
Definition: TICLDumper.cc:243
FreeTrajectoryState outerFreeState(const reco::Track &tk, const MagneticField *field, bool withErr=true)
std::vector< float > stsSC_trackster_barycenter_phi
Definition: TICLDumper.cc:208
std::vector< std::vector< float > > stsCP_trackster_vertices_energy
Definition: TICLDumper.cc:297
std::vector< float > tracksters_merged_barycenter_z
Definition: TICLDumper.cc:343
unsigned int ntracksters_
Definition: TICLDumper.cc:155
const edm::EDGetTokenT< std::vector< reco::Track > > tracks_token_
Definition: TICLDumper.cc:87
std::vector< float > trackster_timeError
Definition: TICLDumper.cc:163
std::vector< std::vector< float > > tracksters_merged_vertices_x
Definition: TICLDumper.cc:356
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeometry_token_
Definition: TICLDumper.cc:108
TTree * trackster_tree_
Definition: TICLDumper.cc:444
edm::ESHandle< MagneticField > bfield_
Definition: TICLDumper.cc:136
~TICLDumper() override
Definition: TICLDumper.cc:815
const edm::EDGetTokenT< std::vector< reco::Muon > > muons_token_
Definition: TICLDumper.cc:105
std::vector< float > stsCP_track_boundaryY
Definition: TICLDumper.cc:283
unsigned int lastLayerEE(bool nose=false) const
Definition: RecHitTools.h:76
std::vector< float > stsCP_boundaryPy
Definition: TICLDumper.cc:280
std::vector< float > simTICLCandidate_time
Definition: TICLDumper.cc:312
std::vector< float > tracksters_merged_barycenter_x
Definition: TICLDumper.cc:341
std::vector< float > stsSC_boundaryEta
Definition: TICLDumper.cc:224
TTree * tree_
Definition: TICLDumper.cc:149
Definition: event.py:1
Definition: Run.h:45
unsigned int ev_event_
Definition: TICLDumper.cc:154
bool saveCLUE3DTracksters_
Definition: TICLDumper.cc:139
std::vector< int > stsSC_trackIdx
Definition: TICLDumper.cc:219
std::vector< float > cluster_position_phi
Definition: TICLDumper.cc:411
std::vector< float > stsCP_track_boundaryEta
Definition: TICLDumper.cc:285
std::vector< float > stsSC_trackster_sigmaPCA2
Definition: TICLDumper.cc:216
const edm::EDGetTokenT< std::vector< SimCluster > > simclusters_token_
Definition: TICLDumper.cc:123
const std::string detector_
Definition: TICLDumper.cc:127
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:381
const edm::EDGetTokenT< ticl::RecoToSimCollectionSimTracksters > MergeRecoToSimSC_token_
Definition: TICLDumper.cc:117
std::vector< float > stsCP_trackster_raw_pt
Definition: TICLDumper.cc:255
std::vector< std::vector< float > > tracksters_merged_vertices_correctedEnergy
Definition: TICLDumper.cc:362
std::vector< float > simTICLCandidate_boundaryPz
Definition: TICLDumper.cc:310