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