CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MuonTrackValidator.h
Go to the documentation of this file.
1 #ifndef MuonTrackValidator_h
2 #define MuonTrackValidator_h
3 
12 
17 
19  public:
22  dirName_ = pset.getParameter<std::string>("dirName");
23  associatormap = pset.getParameter< edm::InputTag >("associatormap");
24  UseAssociators = pset.getParameter< bool >("UseAssociators");
25  tpSelector = TrackingParticleSelector(pset.getParameter<double>("ptMinTP"),
26  pset.getParameter<double>("minRapidityTP"),
27  pset.getParameter<double>("maxRapidityTP"),
28  pset.getParameter<double>("tipTP"),
29  pset.getParameter<double>("lipTP"),
30  pset.getParameter<int>("minHitTP"),
31  pset.getParameter<bool>("signalOnlyTP"),
32  pset.getParameter<bool>("chargedOnlyTP"),
33  pset.getParameter<bool>("stableOnlyTP"),
34  pset.getParameter<std::vector<int> >("pdgIdTP"));
36  pset.getParameter<double>("minRapidityTP"),
37  pset.getParameter<double>("maxRapidityTP"),
38  pset.getParameter<double>("tipTP"),
39  pset.getParameter<double>("lipTP"),
40  pset.getParameter<int>("minHitTP"),
41  pset.getParameter<bool>("chargedOnlyTP"),
42  pset.getParameter<std::vector<int> >("pdgIdTP"));
43 
44  minPhi = pset.getParameter<double>("minPhi");
45  maxPhi = pset.getParameter<double>("maxPhi");
46  nintPhi = pset.getParameter<int>("nintPhi");
47  useGsf = pset.getParameter<bool>("useGsf");
48  BiDirectional_RecoToSim_association = pset.getParameter<bool>("BiDirectional_RecoToSim_association");
49 
50  // dump cfg parameters
51  edm::LogVerbatim("MuonTrackValidator") << "constructing MuonTrackValidator: " << pset.dump();
52 
53  // Declare consumes (also for the base class)
54  bsSrc_Token = consumes<reco::BeamSpot>(bsSrc);
55  tp_effic_Token = consumes<TrackingParticleCollection>(label_tp_effic);
56  tp_fake_Token = consumes<TrackingParticleCollection>(label_tp_fake);
57  for (unsigned int www=0;www<label.size();www++){
59  }
60  simToRecoCollection_Token = consumes<reco::SimToRecoCollection>(associatormap);
61  recoToSimCollection_Token = consumes<reco::RecoToSimCollection>(associatormap);
62 
63  _simHitTpMapTag = mayConsume<SimHitTPAssociationProducer::SimHitTPAssociationList>(pset.getParameter<edm::InputTag>("simHitTpMapTag"));
64 
65  MABH = false;
66  if (!UseAssociators) {
67  // flag MuonAssociatorByHits
68  if (associators[0] == "MuonAssociationByHits") MABH = true;
69  // reset string associators to the map label
70  associators.clear();
71  associators.push_back(associatormap.label());
72  edm::LogVerbatim("MuonTrackValidator") << "--> associators reset to: " <<associators[0];
73  } else {
74  for (auto const& associator :associators) {
75  consumes<reco::TrackToTrackingParticleAssociator>(edm::InputTag(associator));
76  }
77  }
78 
79  // inform on which SimHits will be counted
80  if (usetracker) edm::LogVerbatim("MuonTrackValidator")
81  <<"\n usetracker = TRUE : Tracker SimHits WILL be counted";
82  else edm::LogVerbatim("MuonTrackValidator")
83  <<"\n usetracker = FALSE : Tracker SimHits WILL NOT be counted";
84  if (usemuon) edm::LogVerbatim("MuonTrackValidator")
85  <<" usemuon = TRUE : Muon SimHits WILL be counted";
86  else edm::LogVerbatim("MuonTrackValidator")
87  <<" usemuon = FALSE : Muon SimHits WILL NOT be counted"<<std::endl;
88 
89  // loop over the reco::Track collections to validate: check for inconsistent input settings
90  for (unsigned int www=0;www<label.size();www++) {
91  std::string recoTracksLabel = label[www].label();
92  std::string recoTracksInstance = label[www].instance();
93 
94  // tracks with hits only on tracker
95  if (recoTracksLabel=="generalTracks" ||
96  (recoTracksLabel.find("cutsRecoTracks") != std::string::npos) ||
97  recoTracksLabel=="ctfWithMaterialTracksP5LHCNavigation" ||
98  recoTracksLabel=="hltL3TkTracksFromL2" ||
99  (recoTracksLabel=="hltL3Muons" && recoTracksInstance=="L2Seeded"))
100  {
101  if (usemuon) {
102  edm::LogWarning("MuonTrackValidator")
103  <<"\n*** WARNING : inconsistent input tracksTag = "<<label[www]
104  <<"\n with usemuon == true"<<"\n ---> please change to usemuon == false ";
105  }
106  if (!usetracker) {
107  edm::LogWarning("MuonTrackValidator")
108  <<"\n*** WARNING : inconsistent input tracksTag = "<<label[www]
109  <<"\n with usetracker == false"<<"\n ---> please change to usetracker == true ";
110  }
111  }
112 
113  // tracks with hits only on muon detectors
114  else if (recoTracksLabel=="standAloneMuons" ||
115  recoTracksLabel=="standAloneSETMuons" ||
116  recoTracksLabel=="cosmicMuons" ||
117  recoTracksLabel=="hltL2Muons")
118  {
119  if (usetracker) {
120  edm::LogWarning("MuonTrackValidator")
121  <<"\n*** WARNING : inconsistent input tracksTag = "<<label[www]
122  <<"\n with usetracker == true"<<"\n ---> please change to usetracker == false ";
123  }
124  if (!usemuon) {
125  edm::LogWarning("MuonTrackValidator")
126  <<"\n*** WARNING : inconsistent input tracksTag = "<<label[www]
127  <<"\n with usemuon == false"<<"\n ---> please change to usemuon == true ";
128  }
129  }
130 
131  } // for (unsigned int www=0;www<label.size();www++)
132  }
133 
135  virtual ~MuonTrackValidator(){ }
136 
138  // void beginRun(edm::Run const&, edm::EventSetup const&);
140  void analyze(const edm::Event&, const edm::EventSetup& );
142  void endRun(edm::Run const&, edm::EventSetup const&);
144 
145 private:
147  void getRecoMomentum (const reco::Track& track, double& pt, double& ptError,
148  double& qoverp, double& qoverpError, double& lambda, double& lambdaError,
149  double& phi, double& phiError ) const;
151  void getRecoMomentum (const reco::GsfTrack& gsfTrack, double& pt, double& ptError,
152  double& qoverp, double& qoverpError, double& lambda, double& lambdaError,
153  double& phi, double& phiError) const;
154 
155  private:
161 
163  double minPhi, maxPhi;
164  int nintPhi;
165  bool useGsf;
166  // select tracking particles
167  //(i.e. "denominator" of the efficiency ratio)
170 
171  // flag new validation logic (bidirectional RecoToSim association)
173  // flag MuonAssociatorByHits
174  bool MABH;
175 
176  //1D
177  std::vector<MonitorElement*> h_nchi2, h_nchi2_prob, h_losthits;
178 
179  //2D
180  std::vector<MonitorElement*> chi2_vs_nhits, etares_vs_eta;
181  std::vector<MonitorElement*> h_ptshifteta;
182  std::vector<MonitorElement*> ptres_vs_phi, chi2_vs_phi, nhits_vs_phi, phires_vs_phi;
183 
184  //Profile2D
185  std::vector<MonitorElement*> ptmean_vs_eta_phi, phimean_vs_eta_phi;
186 
187  //assoc chi2
188  std::vector<MonitorElement*> h_assochi2, h_assochi2_prob;
189 
190  //chi2 and # lost hits vs eta: to be used with doProfileX
191  std::vector<MonitorElement*> chi2_vs_eta, nlosthits_vs_eta;
192  std::vector<MonitorElement*> h_chi2meanh, h_losthits_eta;
193  std::vector<MonitorElement*> h_hits_phi;
194  std::vector<MonitorElement*> h_chi2meanhitsh, h_chi2mean_vs_phi;
195 
196  //resolution of track params: to be used with fitslicesytool
199 
200  //pulls of track params vs eta: to be used with fitslicesytool
202  std::vector<MonitorElement*> ptpull_vs_phi, phipull_vs_phi, thetapull_vs_phi;
204  std::vector<MonitorElement*> h_ptpullphi, h_phipullphi, h_thetapullphi;
205 
206 };
207 
208 
209 #endif
std::vector< MonitorElement * > h_dxypulleta
T getParameter(std::string const &) const
std::vector< MonitorElement * > h_ptpulleta
std::vector< MonitorElement * > etares_vs_eta
std::vector< MonitorElement * > ptres_vs_phi
std::vector< edm::EDGetTokenT< edm::View< reco::Track > > > track_Collection_Token
std::vector< MonitorElement * > h_losthits_eta
std::string dump(unsigned int indent=0) const
std::vector< MonitorElement * > h_thetapullphi
std::vector< MonitorElement * > h_ptpullphi
std::vector< MonitorElement * > ptres_vs_eta
std::vector< MonitorElement * > cotThetares_vs_eta
std::vector< MonitorElement * > h_dzpulleta
std::vector< MonitorElement * > h_thetapulleta
std::vector< MonitorElement * > phimean_vs_eta_phi
virtual ~MuonTrackValidator()
Destructor.
std::vector< MonitorElement * > phires_vs_pt
std::vector< MonitorElement * > h_assochi2_prob
std::vector< MonitorElement * > cotThetares_vs_pt
std::vector< MonitorElement * > h_ptshifteta
std::vector< MonitorElement * > dxyres_vs_eta
std::vector< MonitorElement * > dzres_vs_pt
std::vector< MonitorElement * > dzres_vs_eta
void getRecoMomentum(const reco::Track &track, double &pt, double &ptError, double &qoverp, double &qoverpError, double &lambda, double &lambdaError, double &phi, double &phiError) const
retrieval of reconstructed momentum components from reco::Track (== mean values for GSF) ...
std::vector< MonitorElement * > thetapull_vs_phi
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::EDGetTokenT< SimHitTPAssociationProducer::SimHitTPAssociationList > _simHitTpMapTag
std::vector< MonitorElement * > h_assochi2
std::vector< MonitorElement * > chi2_vs_phi
std::vector< MonitorElement * > h_phipulleta
std::vector< MonitorElement * > chi2_vs_nhits
std::vector< MonitorElement * > phires_vs_phi
std::vector< MonitorElement * > h_chi2mean_vs_phi
std::vector< MonitorElement * > h_phipullphi
SingleObjectSelector< TrackingParticleCollection,::TrackingParticleSelector > TrackingParticleSelector
CosmicTrackingParticleSelector cosmictpSelector
std::vector< MonitorElement * > nhits_vs_phi
edm::EDGetTokenT< TrackingParticleCollection > tp_fake_Token
std::vector< MonitorElement * > dxyres_vs_pt
std::vector< edm::InputTag > label
std::vector< MonitorElement * > h_nchi2
std::vector< MonitorElement * > h_nchi2_prob
std::vector< MonitorElement * > ptpull_vs_eta
std::vector< MonitorElement * > dxypull_vs_eta
void analyze(const edm::Event &, const edm::EventSetup &)
Method called before the event loop.
std::vector< MonitorElement * > phipull_vs_eta
std::vector< MonitorElement * > thetapull_vs_eta
MuonTrackValidator(const edm::ParameterSet &pset)
Constructor.
ObjectSelector< CosmicTrackingParticleSelector > CosmicTrackingParticleSelector
std::vector< MonitorElement * > ptmean_vs_eta_phi
edm::EDGetTokenT< reco::RecoToSimCollection > recoToSimCollection_Token
edm::EDGetTokenT< reco::SimToRecoCollection > simToRecoCollection_Token
std::vector< MonitorElement * > chi2_vs_eta
TrackingParticleSelector tpSelector
std::vector< MonitorElement * > h_losthits
std::vector< MonitorElement * > ptpull_vs_phi
edm::InputTag associatormap
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &)
std::vector< MonitorElement * > ptres_vs_pt
std::vector< MonitorElement * > h_chi2meanhitsh
std::vector< MonitorElement * > phires_vs_eta
std::vector< MonitorElement * > h_hits_phi
std::vector< std::string > associators
std::vector< MonitorElement * > nlosthits_vs_eta
std::vector< MonitorElement * > dzpull_vs_eta
edm::EDGetTokenT< TrackingParticleCollection > tp_effic_Token
Definition: Run.h:41
std::vector< MonitorElement * > h_chi2meanh
void endRun(edm::Run const &, edm::EventSetup const &)
Method called at the end of the event loop.
std::vector< MonitorElement * > phipull_vs_phi
edm::EDGetTokenT< reco::BeamSpot > bsSrc_Token
Definition: DDAxes.h:10