CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
QualityCutsAnalyzer.cc
Go to the documentation of this file.
1 #include <algorithm>
2 #include <cctype>
3 #include <iomanip>
4 #include <set>
5 #include <sstream>
6 #include <vector>
7 
8 #include "TFile.h"
9 #include "TH1F.h"
10 
11 #include "HepPDT/ParticleID.hh"
12 
13 // user include files
17 
23 
25 
27 
29 
31 
32 //
33 // class decleration
34 //
35 
37 {
38 
39 public:
40 
41  explicit QualityCutsAnalyzer(const edm::ParameterSet&);
43 
44 private:
45 
46  virtual void beginJob() override ;
47  virtual void analyze(const edm::Event&, const edm::EventSetup&) override;
48  virtual void endJob() override;
49 
50  // Member data
51 
52  typedef std::vector<int> vint;
53  typedef std::vector<std::string> vstring;
54 
58 
60 
63 
66 
67  void
72  );
73 
74  // Histograms for optimization
75 
77  {
78  double sdl; // Signed decay length
79  double dta; // Distance to jet axis
80  double tip; // Transverse impact parameter
81  double lip; // Longitudinal impact parameter
82  double ips; // Impact parameter significance.
83  double pt; // Transverse momentum
84  double chi2; // Chi^2
85  std::size_t hits; // Number of hits
86  std::size_t pixelhits; // Number of hits
87 
88  histogram_element_t(double d, double a, double t, double l, double i, double p, double c, std::size_t h, std::size_t x)
89  {
90  sdl = d;
91  dta = a;
92  tip = t;
93  lip = l;
94  ips = i;
95  pt = p;
96  chi2 = c;
97  hits = h;
98  pixelhits = x;
99  }
100 
102  {
103  sdl = orig.sdl;
104  dta = orig.dta;
105  tip = orig.tip;
106  lip = orig.lip;
107  ips = orig.ips;
108  pt = orig.pt;
109  chi2 = orig.chi2;
110  hits = orig.hits;
111  pixelhits = orig.pixelhits;
112  }
113  };
114 
115  typedef std::vector<std::vector<histogram_element_t> > histogram_data_t;
117 
119  {
120 
121  TH1F* sdl;
122  TH1F* dta;
123  TH1F* tip;
124  TH1F* lip;
125  TH1F* ips;
126  TH1F* pixelhits;
127  TH1F* pt_1gev;
128  TH1F* chi2;
129  TH1F* hits;
130 
131  public:
132 
134  {
136  name = std::string("hits_") + particleType;
137  title = std::string("Hit distribution for ") + particleType;
138  hits = new TH1F(name.c_str(), title.c_str(), 19, -0.5, 18.5);
139 
140  name = std::string("chi2_") + particleType;
141  title = std::string("Chi2 distribution for ") + particleType;
142  chi2 = new TH1F(name.c_str(), title.c_str(), 100, 0., 30.);
143 
144  name = std::string("pixelhits_") + particleType;
145  title = std::string("Pixel hits distribution for ") + particleType;
146  pixelhits = new TH1F(name.c_str(), title.c_str(), 21, -0.5, 20.5);
147 
148  name = std::string("pt_1Gev_") + particleType;
149  title = std::string("Pt distribution close 1Gev for ") + particleType;
150  pt_1gev = new TH1F(name.c_str(), title.c_str(), 100, 0., 2.);
151 
152  name = std::string("tip_") + particleType;
153  title = std::string("Transverse impact parameter distribution for ") + particleType;
154  tip = new TH1F(name.c_str(), title.c_str(), 100, -0.3, 0.3);
155 
156  name = std::string("lip_") + particleType;
157  title = std::string("Longitudinal impact parameter distribution for ") + particleType;
158  lip = new TH1F(name.c_str(), title.c_str(), 100, -1., 1.);
159 
160  name = std::string("ips_") + particleType;
161  title = std::string("IPS distribution for ") + particleType;
162  ips = new TH1F(name.c_str(), title.c_str(), 100, -25.0, 25.0);
163 
164  name = std::string("sdl_") + particleType;
165  title = std::string("Decay length distribution for ") + particleType;
166  sdl = new TH1F(name.c_str(), title.c_str(), 100, -5., 5.);
167 
168  name = std::string("dta_") + particleType;
169  title = std::string("Distance to jet distribution for ") + particleType;
170  dta = new TH1F(name.c_str(), title.c_str(), 100, 0.0, 0.2);
171  }
172 
174  {
175  delete hits;
176  delete chi2;
177  delete pixelhits;
178  delete pt_1gev;
179  delete tip;
180  delete lip;
181  delete ips;
182  delete sdl;
183  delete dta;
184  }
185 
187  {
188  hits->Fill(data.hits);
189  chi2->Fill(data.chi2);
190  pixelhits->Fill(data.pt);
191  pt_1gev->Fill(data.pt);
192  ips->Fill(data.ips);
193  tip->Fill(data.tip);
194  lip->Fill(data.lip);
195  sdl->Fill(data.sdl);
196  dta->Fill(data.dta);
197  }
198 
199  void Write()
200  {
201  hits->Write();
202  chi2->Write();
203  pixelhits->Write();
204  pt_1gev->Write();
205  ips->Write();
206  tip->Write();
207  lip->Write();
208  sdl->Write();
209  dta->Write();
210  }
211  };
212 
213  // Track classification.
215 
216 };
217 
218 
219 //
220 // constructors and destructor
221 //
222 QualityCutsAnalyzer::QualityCutsAnalyzer(const edm::ParameterSet& config) : classifier_(config,consumesCollector())
223 {
224  trackProducer_ = config.getUntrackedParameter<edm::InputTag> ( "trackProducer" );
225  consumes<edm::View<reco::Track>>(trackProducer_);
226  primaryVertexProducer_ = config.getUntrackedParameter<edm::InputTag> ( "primaryVertexProducer" );
227  consumes<reco::VertexCollection>(primaryVertexProducer_);
228  jetTracksAssociation_ = config.getUntrackedParameter<edm::InputTag> ( "jetTracksAssociation" );
229  consumes<reco::JetTracksAssociationCollection>(jetTracksAssociation_);
230 
231  rootFile_ = config.getUntrackedParameter<std::string> ( "rootFile" );
232 
233  minimumNumberOfHits_ = config.getUntrackedParameter<int> ( "minimumNumberOfHits" );
234  minimumNumberOfPixelHits_ = config.getUntrackedParameter<int> ( "minimumNumberOfPixelHits" );
235  minimumTransverseMomentum_ = config.getUntrackedParameter<double> ( "minimumTransverseMomentum" );
236  maximumChiSquared_ = config.getUntrackedParameter<double> ( "maximumChiSquared" );
237 
238  std::string trackQualityType = config.getUntrackedParameter<std::string>("trackQualityClass"); //used
239  trackQuality_ = reco::TrackBase::qualityByName(trackQualityType);
240  useAllQualities_ = false;
241 
242  std::transform(trackQualityType.begin(), trackQualityType.end(), trackQualityType.begin(), (int(*)(int)) std::tolower);
243  if (trackQualityType == "any")
244  {
245  std::cout << "Using any" << std::endl;
246  useAllQualities_ = true;
247  }
248 }
249 
251 
252 //
253 // member functions
254 //
255 
256 // ------------ method called to for each event ------------
257 void
259 {
260  // Track collection
262  event.getByLabel(trackProducer_,trackCollection);
263  // Primary vertex
264  edm::Handle<reco::VertexCollection> primaryVertexCollection;
265  event.getByLabel(primaryVertexProducer_, primaryVertexCollection);
266  // Jet to tracks associator
268  event.getByLabel(jetTracksAssociation_, jetTracks);
269  // Trasient track builder
271  setup.get<TransientTrackRecord>().get("TransientTrackBuilder", TTbuilder);
272 
273  // Setting up event information for the track categories.
274  classifier_.newEvent(event, setup);
275 
277  TTbuilder,
278  primaryVertexCollection,
279  jetTracks
280  );
281 }
282 
283 
284 // ------------ method called once each job just before starting event loop ------------
285 void
287 {
288  histogram_data_.resize(6);
289 }
290 
291 
292 // ------------ method called once each job just after ending the event loop ------------
293 void
295 {
296  TFile file(rootFile_.c_str(), "RECREATE");
297  file.cd();
298 
299  // saving the histograms
300  for (std::size_t i=0; i<6; i++)
301  {
302  std::string particle;
303  if (i == 0)
304  particle = std::string("B_tracks");
305  else if (i == 1)
306  particle = std::string("C_tracks");
307  else if (i == 2)
308  particle = std::string("nonB_tracks");
309  else if (i == 3)
310  particle = std::string("displaced_tracks");
311  else if (i == 4)
312  particle = std::string("bad_tracks");
313  else
314  particle = std::string("fake_tracks");
315 
316  histogram_t histogram(particle);
317 
318  for (std::size_t j=0; j<histogram_data_[i].size(); j++)
319  histogram.Fill(histogram_data_[i][j]);
320 
321  histogram.Write();
322  }
323 
324  file.Flush();
325 }
326 
327 
328 void
330  const edm::ESHandle<TransientTrackBuilder> & TTbuilder,
331  const edm::Handle<reco::VertexCollection> & primaryVertexProducer_,
332  const edm::Handle<reco::JetTracksAssociationCollection> & jetTracksAssociation
333 )
334 {
335  const TransientTrackBuilder * bproduct = TTbuilder.product();
336 
337  // getting the primary vertex
338  // use first pv of the collection
340 
341  if (primaryVertexProducer_->size() != 0)
342  {
344  std::vector<reco::Vertex> sortedList = pvs.sortedList(*(primaryVertexProducer_.product()));
345  pv = (sortedList.front());
346  }
347  else
348  { // create a dummy PV
349  // cout << "NO PV FOUND" << endl;
351  e(0,0)=0.0015*0.0015;
352  e(1,1)=0.0015*0.0015;
353  e(2,2)=15.*15.;
354  reco::Vertex::Point p(0,0,0);
355  pv = reco::Vertex(p,e,1,1,1);
356  }
357 
358  reco::JetTracksAssociationCollection::const_iterator it = jetTracksAssociation->begin();
359 
360  int i=0;
361 
362  for (; it != jetTracksAssociation->end(); it++, i++)
363  {
364  // get jetTracks object
365  reco::JetTracksAssociationRef jetTracks(jetTracksAssociation, i);
366 
367  double pvZ = pv.z();
368  GlobalVector direction(jetTracks->first->px(), jetTracks->first->py(), jetTracks->first->pz());
369 
370  // get the tracks associated to the jet
371  reco::TrackRefVector tracks = jetTracks->second;
372  for (std::size_t index = 0; index < tracks.size(); index++)
373  {
374  edm::RefToBase<reco::Track> track(tracks[index]);
375 
376  double pt = tracks[index]->pt();
377  double chi2 = tracks[index]->normalizedChi2();
378  int hits = tracks[index]->hitPattern().numberOfValidHits();
379  int pixelHits = tracks[index]->hitPattern().numberOfValidPixelHits();
380 
381  if (
382  hits < minimumNumberOfHits_ ||
383  pixelHits < minimumNumberOfPixelHits_ ||
385  chi2 > maximumChiSquared_ ||
386  (!useAllQualities_ && !tracks[index]->quality(trackQuality_))
387  ) continue;
388 
389  const reco::TransientTrack transientTrack = bproduct->build(&(*tracks[index]));
390  double dta = - IPTools::jetTrackDistance(transientTrack, direction, pv).second.value();
391  double sdl = IPTools::signedDecayLength3D(transientTrack, direction, pv).second.value();
392  double ips = IPTools::signedImpactParameter3D(transientTrack, direction, pv).second.value();
393  double d0 = IPTools::signedTransverseImpactParameter(transientTrack, direction, pv).second.value();
394  double dz = tracks[index]->dz() - pvZ;
395 
396  // Classify the reco track;
398 
399  // Check for the different categories
401  histogram_data_[5].push_back(histogram_element_t(sdl, dta, d0, dz, ips, pt, chi2, hits, pixelHits));
403  histogram_data_[0].push_back(histogram_element_t(sdl, dta, d0, dz, ips, pt, chi2, hits, pixelHits));
404  else if ( classifier_.is(TrackClassifier::Bad) )
405  histogram_data_[4].push_back(histogram_element_t(sdl, dta, d0, dz, ips, pt, chi2, hits, pixelHits));
406  else if (
409  )
410  histogram_data_[3].push_back(histogram_element_t(sdl, dta, d0, dz, ips, pt, chi2, hits, pixelHits));
412  histogram_data_[1].push_back(histogram_element_t(sdl, dta, d0, dz, ips, pt, chi2, hits, pixelHits));
413  else
414  histogram_data_[2].push_back(histogram_element_t(sdl, dta, d0, dz, ips, pt, chi2, hits, pixelHits));
415 
416  }
417  }
418 }
419 
T getUntrackedParameter(std::string const &, T const &) const
tuple t
Definition: tree.py:139
int i
Definition: DBlmapReader.cc:9
histogram_element_t(double d, double a, double t, double l, double i, double p, double c, std::size_t h, std::size_t x)
std::vector< std::vector< histogram_element_t > > histogram_data_t
transient_vector_type::const_iterator const_iterator
TrackQuality
track quality
Definition: TrackBase.h:149
std::vector< int > vint
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
reco::TrackBase::TrackQuality trackQuality_
std::pair< bool, Measurement1D > signedTransverseImpactParameter(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:50
void LoopOverJetTracksAssociation(const edm::ESHandle< TransientTrackBuilder > &, const edm::Handle< reco::VertexCollection > &, const edm::Handle< reco::JetTracksAssociationCollection > &)
std::pair< bool, Measurement1D > signedDecayLength3D(const TrajectoryStateOnSurface &state, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:94
reco::TransientTrack build(const reco::Track *p) const
std::vector< reco::Vertex > sortedList(const reco::VertexCollection &primaryVertex) const
virtual void endJob() override
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:43
std::pair< bool, Measurement1D > signedImpactParameter3D(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:71
TrackClassifier classifier_
std::pair< double, Measurement1D > jetTrackDistance(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:200
histogram_element_t(const histogram_element_t &orig)
tuple d
Definition: ztail.py:151
TrackClassifier const & evaluate(reco::TrackBaseRef const &)
Classify the RecoTrack in categories.
histogram_t(const std::string &particleType)
int j
Definition: DBlmapReader.cc:9
double z() const
y coordinate
Definition: Vertex.h:112
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
void newEvent(edm::Event const &, edm::EventSetup const &)
Pre-process event information (for accessing reconstraction information)
edm::InputTag trackProducer_
bool is(Category category) const
Returns track flag for a given category.
std::vector< std::string > vstring
Get track history and classify it in function of their .
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:123
T const * product() const
Definition: Handle.h:81
tuple tracks
Definition: testEve_cfg.py:39
const T & get() const
Definition: EventSetup.h:56
T const * product() const
Definition: ESHandle.h:86
void Fill(const histogram_element_t &data)
virtual void beginJob() override
QualityCutsAnalyzer(const edm::ParameterSet &)
double a
Definition: hdecay.h:121
edm::InputTag primaryVertexProducer_
tuple cout
Definition: gather_cfg.py:121
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
histogram_data_t histogram_data_
edm::InputTag jetTracksAssociation_