test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MultiTrackValidatorGenPs.cc
Go to the documentation of this file.
3 
6 
23 
27 
28 #include "TMath.h"
29 #include <TF1.h>
30 
31 //#include <iostream>
32 
33 using namespace std;
34 using namespace edm;
35 
36 static const std::string kTrackAssociatorByChi2("trackAssociatorByChi2");
37 
39 
40  gpSelector = GenParticleCustomSelector(pset.getParameter<double>("ptMinGP"),
41  pset.getParameter<double>("minRapidityGP"),
42  pset.getParameter<double>("maxRapidityGP"),
43  pset.getParameter<double>("tipGP"),
44  pset.getParameter<double>("lipGP"),
45  pset.getParameter<bool>("chargedOnlyGP"),
46  pset.getParameter<int>("statusGP"),
47  pset.getParameter<std::vector<int> >("pdgIdGP"));
48 
49  if(UseAssociators) {
50  for(auto const& src: associators) {
51  if( src.label() == kTrackAssociatorByChi2) {
52  label_gen_associator = consumes<reco::TrackToGenParticleAssociator>(src);
53  break;
54  }
55  }
56  }
57  else {
58  for (auto const& src: associators) {
59  associatormapGtR = consumes<reco::GenToRecoCollection>(src);
60  associatormapRtG = consumes<reco::RecoToGenCollection>(src);
61  break;
62  }
63  }
64 }
65 
67 
68 
70  using namespace reco;
71 
72  edm::LogInfo("TrackValidator") << "\n====================================================" << "\n"
73  << "Analyzing new event" << "\n"
74  << "====================================================\n" << "\n";
75 
76  edm::ESHandle<ParametersDefinerForTP> parametersDefinerTP;
77  setup.get<TrackAssociatorRecord>().get(parametersDefiner,parametersDefinerTP);
78 
79  edm::Handle<GenParticleCollection> TPCollectionHeff ;
80  event.getByToken(label_tp_effic,TPCollectionHeff);
81  const GenParticleCollection tPCeff = *(TPCollectionHeff.product());
82 
83  edm::Handle<GenParticleCollection> TPCollectionHfake ;
84  event.getByToken(label_tp_fake,TPCollectionHfake);
85  const GenParticleCollection tPCfake = *(TPCollectionHfake.product());
86 
87  //if (tPCeff.size()==0) {edm::LogInfo("TrackValidator")
88  //<< "TP Collection for efficiency studies has size = 0! Skipping Event." ; return;}
89  //if (tPCfake.size()==0) {edm::LogInfo("TrackValidator")
90  //<< "TP Collection for fake rate studies has size = 0! Skipping Event." ; return;}
91 
92  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
93  event.getByToken(bsSrc,recoBeamSpotHandle);
94  reco::BeamSpot bs = *recoBeamSpotHandle;
95 
97  event.getByToken(label_pileupinfo,puinfoH);
98  PileupSummaryInfo puinfo;
99 
100  for (unsigned int puinfo_ite=0;puinfo_ite<(*puinfoH).size();++puinfo_ite){
101  if ((*puinfoH)[puinfo_ite].getBunchCrossing()==0){
102  puinfo=(*puinfoH)[puinfo_ite];
103  break;
104  }
105  }
106 
107  const reco::TrackToGenParticleAssociator* trackGenAssociator =nullptr;
108  if(UseAssociators) {
110  return;
111  }
112  else {
114  event.getByToken(label_gen_associator,trackGenAssociatorH);
115  trackGenAssociator = trackGenAssociatorH.product();
116  }
117  }
118  else if(associatormapGtR.isUninitialized()) {
119  return;
120  }
121 
122  // dE/dx
123  // at some point this could be generalized, with a vector of tags and a corresponding vector of Handles
124  // I'm writing the interface such to take vectors of ValueMaps
125  std::vector<const edm::ValueMap<reco::DeDxData> *> v_dEdx;
126  if(dodEdxPlots_) {
129  event.getByToken(m_dEdx1Tag, dEdx1Handle);
130  event.getByToken(m_dEdx2Tag, dEdx2Handle);
131  v_dEdx.push_back(dEdx1Handle.product());
132  v_dEdx.push_back(dEdx2Handle.product());
133  }
134 
135  int w=0; //counter counting the number of sets of histograms
136  for (unsigned int www=0;www<label.size();www++){
137  //
138  //get collections from the event
139  //
141  if(!event.getByToken(labelToken[www], trackCollection)&&ignoremissingtkcollection_)continue;
142  //if (trackCollection->size()==0)
143  //edm::LogInfo("TrackValidator") << "TrackCollection size = 0!" ;
144  //continue;
145  //}
146  reco::RecoToGenCollection recGenColl;
147  reco::GenToRecoCollection genRecColl;
148  //associate tracks
149  if(UseAssociators){
150  edm::LogVerbatim("TrackValidator") << "Analyzing "
151  << label[www].process()<<":"
152  << label[www].label()<<":"
153  << label[www].instance()<<" with "
154  << kTrackAssociatorByChi2 <<"\n";
155 
156  LogTrace("TrackValidator") << "Calling associateRecoToGen method" << "\n";
157  recGenColl=trackGenAssociator->associateRecoToGen(trackCollection,
158  TPCollectionHfake);
159  LogTrace("TrackValidator") << "Calling associateGenToReco method" << "\n";
160  genRecColl=trackGenAssociator->associateGenToReco(trackCollection,
161  TPCollectionHeff);
162  }
163  else{
164  edm::LogVerbatim("TrackValidator") << "Analyzing "
165  << label[www].process()<<":"
166  << label[www].label()<<":"
167  << label[www].instance()<<" with "
168  << associators[0] << "\n";
169 
170  Handle<reco::GenToRecoCollection > gentorecoCollectionH;
171  event.getByToken(associatormapGtR,gentorecoCollectionH);
172  genRecColl= *(gentorecoCollectionH.product());
173 
174  Handle<reco::RecoToGenCollection > recotogenCollectionH;
175  event.getByToken(associatormapRtG,recotogenCollectionH);
176  recGenColl= *(recotogenCollectionH.product());
177  }
178 
179 
180 
181  // ########################################################
182  // fill simulation histograms (LOOP OVER TRACKINGPARTICLES)
183  // ########################################################
184 
185  //compute number of tracks per eta interval
186  //
187  edm::LogVerbatim("TrackValidator") << "\n# of GenParticles: " << tPCeff.size() << "\n";
188  int ats(0); //This counter counts the number of simTracks that are "associated" to recoTracks
189  int st(0); //This counter counts the number of simulated tracks passing the MTV selection (i.e. tpSelector(tp) )
190  unsigned sts(0); //This counter counts the number of simTracks surviving the bunchcrossing cut
191  unsigned asts(0); //This counter counts the number of simTracks that are "associated" to recoTracks surviving the bunchcrossing cut
192  for (GenParticleCollection::size_type i=0; i<tPCeff.size(); i++){ //loop over TPs collection for tracking efficiency
193  GenParticleRef tpr(TPCollectionHeff, i);
194  GenParticle* tp=const_cast<GenParticle*>(tpr.get());
195  TrackingParticle::Vector momentumTP;
196  TrackingParticle::Point vertexTP;
197  double dxyGen(0);
198  double dzGen(0);
199 
200 
201  //---------- THIS PART HAS TO BE CLEANED UP. THE PARAMETER DEFINER WAS NOT MEANT TO BE USED IN THIS WAY ----------
202  //If the GenParticle is collison like, get the momentum and vertex at production state
204  {
205  //fixme this one shold be implemented
206  if(! gpSelector(*tp)) continue;
207  momentumTP = tp->momentum();
208  vertexTP = tp->vertex();
209  //Calcualte the impact parameters w.r.t. PCA
210  TrackingParticle::Vector momentum = parametersDefinerTP->momentum(event,setup,*tp);
211  TrackingParticle::Point vertex = parametersDefinerTP->vertex(event,setup,*tp);
212  dxyGen = (-vertex.x()*sin(momentum.phi())+vertex.y()*cos(momentum.phi()));
213  dzGen = vertex.z() - (vertex.x()*momentum.x()+vertex.y()*momentum.y())/sqrt(momentum.perp2())
214  * momentum.z()/sqrt(momentum.perp2());
215  }
216  //If the GenParticle is comics, get the momentum and vertex at PCA
217  else
218  {
219  //if(! cosmictpSelector(*tp,&bs,event,setup)) continue;
220  momentumTP = parametersDefinerTP->momentum(event,setup,*tp);
221  vertexTP = parametersDefinerTP->vertex(event,setup,*tp);
222  dxyGen = (-vertexTP.x()*sin(momentumTP.phi())+vertexTP.y()*cos(momentumTP.phi()));
223  dzGen = vertexTP.z() - (vertexTP.x()*momentumTP.x()+vertexTP.y()*momentumTP.y())/sqrt(momentumTP.perp2())
224  * momentumTP.z()/sqrt(momentumTP.perp2());
225  }
226  //---------- THE PART ABOVE HAS TO BE CLEANED UP. THE PARAMETER DEFINER WAS NOT MEANT TO BE USED IN THIS WAY ----------
227 
228  st++; //This counter counts the number of simulated tracks passing the MTV selection (i.e. tpSelector(tp) )
229 
230  // in the coming lines, histos are filled using as input
231  // - momentumTP
232  // - vertexTP
233  // - dxyGen
234  // - dzGen
235 
236  if(doSimPlots_ && w == 0) {
237  histoProducerAlgo_->fill_generic_simTrack_histos(momentumTP,vertexTP, tp->collisionId());//fixme: check meaning of collisionId
238  }
239  if(!doSimTrackPlots_)
240  continue;
241 
242 
243  // ##############################################
244  // fill RecoAssociated GenTracks' histograms
245  // ##############################################
246  // bool isRecoMatched(false); // UNUSED
247  const reco::Track* matchedTrackPointer=0;
248  std::vector<std::pair<RefToBase<Track>, double> > rt;
249  if(genRecColl.find(tpr) != genRecColl.end()){
250  rt = (std::vector<std::pair<RefToBase<Track>, double> >) genRecColl[tpr];
251  if (rt.size()!=0) {
252  ats++; //This counter counts the number of simTracks that have a recoTrack associated
253  // isRecoMatched = true; // UNUSED
254  matchedTrackPointer = rt.begin()->first.get();
255  edm::LogVerbatim("TrackValidator") << "GenParticle #" << st
256  << " with pt=" << sqrt(momentumTP.perp2())
257  << " associated with quality:" << rt.begin()->second <<"\n";
258  }
259  }else{
260  edm::LogVerbatim("TrackValidator")
261  << "GenParticle #" << st
262  << " with pt,eta,phi: "
263  << sqrt(momentumTP.perp2()) << " , "
264  << momentumTP.eta() << " , "
265  << momentumTP.phi() << " , "
266  << " NOT associated to any reco::Track" << "\n";
267  }
268 
269 
270 
271  int nSimHits = 0;
272  histoProducerAlgo_->fill_recoAssociated_simTrack_histos(w,*tp,momentumTP,vertexTP,dxyGen,dzGen,nSimHits,matchedTrackPointer,puinfo.getPU_NumInteractions());
273 
274  sts++;
275  if (matchedTrackPointer) asts++;
276 
277 
278 
279 
280  } // End for (GenParticleCollection::size_type i=0; i<tPCeff.size(); i++){
281 
282  if(doSimPlots_ && w == 0) {
283  histoProducerAlgo_->fill_simTrackBased_histos(st);
284  }
285 
286 
287  // ##############################################
288  // fill recoTracks histograms (LOOP OVER TRACKS)
289  // ##############################################
290  if(!doRecoTrackPlots_)
291  continue;
292  edm::LogVerbatim("TrackValidator") << "\n# of reco::Tracks with "
293  << label[www].process()<<":"
294  << label[www].label()<<":"
295  << label[www].instance()
296  << ": " << trackCollection->size() << "\n";
297 
298  //int sat(0); //This counter counts the number of recoTracks that are associated to GenTracks from Signal only
299  int at(0); //This counter counts the number of recoTracks that are associated to GenTracks
300  int rT(0); //This counter counts the number of recoTracks in general
301 
302  for(View<Track>::size_type i=0; i<trackCollection->size(); ++i){
303 
304  RefToBase<Track> track(trackCollection, i);
305  rT++;
306 
307  bool isSigGenMatched(false);
308  bool isGenMatched(false);
309  bool isChargeMatched(true);
310  int numAssocRecoTracks = 0;
311  int nSimHits = 0;
312  double sharedFraction = 0.;
313  std::vector<std::pair<GenParticleRef, double> > tp;
314  if(recGenColl.find(track) != recGenColl.end()){
315  tp = recGenColl[track];
316  if (tp.size()!=0) {
317  /*
318  std::vector<PSimHit> simhits=tp[0].first->trackPSimHit(DetId::Tracker);
319  nSimHits = simhits.end()-simhits.begin();
320  */
321  sharedFraction = tp[0].second;
322  isGenMatched = true;
323  if (tp[0].first->charge() != track->charge()) isChargeMatched = false;
324  if(genRecColl.find(tp[0].first) != genRecColl.end()) numAssocRecoTracks = genRecColl[tp[0].first].size();
325  //std::cout << numAssocRecoTracks << std::endl;
326  at++;
327  for (unsigned int tp_ite=0;tp_ite<tp.size();++tp_ite){
328  GenParticle trackpart = *(tp[tp_ite].first);
329  /*
330  if ((trackpart.eventId().event() == 0) && (trackpart.eventId().bunchCrossing() == 0)){
331  isSigGenMatched = true;
332  sat++;
333  break;
334  }
335  */
336  }
337  edm::LogVerbatim("TrackValidator") << "reco::Track #" << rT << " with pt=" << track->pt()
338  << " associated with quality:" << tp.begin()->second <<"\n";
339  }
340  } else {
341  edm::LogVerbatim("TrackValidator") << "reco::Track #" << rT << " with pt=" << track->pt()
342  << " NOT associated to any GenParticle" << "\n";
343  }
344 
345 
346  double dR=0;//fixme: plots vs dR not implemented for now
347  histoProducerAlgo_->fill_generic_recoTrack_histos(w,*track,bs.position(), nullptr, isGenMatched,isSigGenMatched, isChargeMatched, numAssocRecoTracks, puinfo.getPU_NumInteractions(), nSimHits, sharedFraction,dR);
348 
349  // dE/dx
350  if (dodEdxPlots_) histoProducerAlgo_->fill_dedx_recoTrack_histos(w,track, v_dEdx);
351 
352 
353  //Fill other histos
354  //try{ //Is this really necessary ????
355 
356  if (tp.size()==0) continue;
357 
358  histoProducerAlgo_->fill_simAssociated_recoTrack_histos(w,*track);
359 
360  GenParticleRef tpr = tp.begin()->first;
361 
362  /* TO BE FIXED LATER
363  if (associators[ww]=="TrackAssociatorByChi2"){
364  //association chi2
365  double assocChi2 = -tp.begin()->second;//in association map is stored -chi2
366  h_assochi2[www]->Fill(assocChi2);
367  h_assochi2_prob[www]->Fill(TMath::Prob((assocChi2)*5,5));
368  }
369  else if (associators[ww]=="quickTrackAssociatorByHits"){
370  double fraction = tp.begin()->second;
371  h_assocFraction[www]->Fill(fraction);
372  h_assocSharedHit[www]->Fill(fraction*track->numberOfValidHits());
373  }
374  */
375 
376 
377  //Get tracking particle parameters at point of closest approach to the beamline
378  TrackingParticle::Vector momentumTP = parametersDefinerTP->momentum(event,setup,*(tpr.get()));
379  TrackingParticle::Point vertexTP = parametersDefinerTP->vertex(event,setup,*(tpr.get()));
380  int chargeTP = tpr->charge();
381 
382  histoProducerAlgo_->fill_ResoAndPull_recoTrack_histos(w,momentumTP,vertexTP,chargeTP,
383  *track,bs.position());
384 
385 
386  //TO BE FIXED
387  //std::vector<PSimHit> simhits=tpr.get()->trackPSimHit(DetId::Tracker);
388  //nrecHit_vs_nsimHit_rec2sim[w]->Fill(track->numberOfValidHits(), (int)(simhits.end()-simhits.begin() ));
389 
390  /*
391  } // End of try{
392  catch (cms::Exception e){
393  LogTrace("TrackValidator") << "exception found: " << e.what() << "\n";
394  }
395  */
396 
397  } // End of for(View<Track>::size_type i=0; i<trackCollection->size(); ++i){
398 
399  histoProducerAlgo_->fill_trackBased_histos(w,at,rT,st);
400 
401  edm::LogVerbatim("TrackValidator") << "Total Simulated: " << st << "\n"
402  << "Total Associated (genToReco): " << ats << "\n"
403  << "Total Reconstructed: " << rT << "\n"
404  << "Total Associated (recoToGen): " << at << "\n"
405  << "Total Fakes: " << rT-at << "\n";
406 
407  w++;
408  } // End of for (unsigned int www=0;www<label.size();www++){
409 }
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
T getParameter(std::string const &) const
unsigned int size_type
Definition: View.h:85
edm::EDGetTokenT< edm::ValueMap< reco::DeDxData > > m_dEdx1Tag
edm::EDGetTokenT< reco::GenToRecoCollection > associatormapGtR
int i
Definition: DBlmapReader.cc:9
std::vector< edm::InputTag > associators
static const std::string kTrackAssociatorByChi2("trackAssociatorByChi2")
const double w
Definition: UKUtility.cc:23
edm::EDGetTokenT< std::vector< PileupSummaryInfo > > label_pileupinfo
const_iterator end() const
last iterator over the map (read only)
edm::EDGetTokenT< TrackingParticleCollection > label_tp_effic
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
virtual const Point & vertex() const
vertex position (overwritten by PF...)
void analyze(const edm::Event &, const edm::EventSetup &)
Method called once per event.
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
const_iterator find(const key_type &k) const
find element with specified reference key
int collisionId() const
Definition: GenParticle.h:39
uint16_t size_type
MultiTrackValidatorGenPs(const edm::ParameterSet &pset)
Constructor.
edm::EDGetTokenT< reco::BeamSpot > bsSrc
math::XYZPointD Point
point in the space
edm::EDGetTokenT< reco::TrackToGenParticleAssociator > label_gen_associator
std::vector< edm::InputTag > label
GenParticleCustomSelector gpSelector
T sqrt(T t)
Definition: SSEVec.h:18
double pt() const
track transverse momentum
Definition: TrackBase.h:616
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:244
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
#define LogTrace(id)
std::unique_ptr< MTVHistoProducerAlgoForTracker > histoProducerAlgo_
virtual Vector momentum() const final
spatial momentum vector
T const * product() const
Definition: Handle.h:81
const T & get() const
Definition: EventSetup.h:56
bool isUninitialized() const
Definition: EDGetToken.h:73
int charge() const
track electric charge
Definition: TrackBase.h:562
SingleObjectSelector< GenParticleCollection,::GenParticleCustomSelector > GenParticleCustomSelector
std::vector< edm::EDGetTokenT< edm::View< reco::Track > > > labelToken
math::XYZVectorD Vector
point in the space
edm::EDGetTokenT< TrackingParticleCollection > label_tp_fake
edm::EDGetTokenT< reco::RecoToGenCollection > associatormapRtG
tuple size
Write out results.
virtual ~MultiTrackValidatorGenPs()
Destructor.
edm::EDGetTokenT< edm::ValueMap< reco::DeDxData > > m_dEdx2Tag
const bool parametersDefinerIsCosmic_