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