ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/Jeng/PVStudy/plugins/PVEffAnalyzer.cc
Revision: 1.8
Committed: Fri Jun 4 17:46:09 2010 UTC (14 years, 11 months ago) by yygao
Content type: text/plain
Branch: MAIN
CVS Tags: CMSSW_3_6_1_patch3_May27thReReco, HEAD
Changes since 1.7: +3 -3 lines
Error occurred while calculating annotation data.
Log Message:
fix bug and switch ClusterizerInZ to GapClusterizer

File Contents

# Content
1
2 // -*- C++ -*-
3 //
4 // Package: PVEffAnalyzer
5 // Class: PVEffAnalyzer
6 //
7 /**\class PVEffAnalyzer PVEffAnalyzer.cc UserCode/PVEffAnalyzer/plugins/PVEffAnalyzer.cc
8
9 Description: <one line class summary>
10
11 Implementation:
12 <Notes on implementation>
13 */
14 //
15 // Original Author: "Geng-yuan Jeng/UC Riverside"
16 // "Yanyan Gao/Fermilab ygao@fnal.gov"
17 // Created: Thu Aug 20 11:55:40 CDT 2009
18 // $Id: PVEffAnalyzer.cc,v 1.7 2010/05/13 09:28:30 yygao Exp $
19 //
20 //
21
22
23 // system include files
24 #include <memory>
25 #include <string>
26 #include <vector>
27 #include <iostream>
28 #include <sstream>
29
30 // user include files
31 #include "FWCore/Framework/interface/Frameworkfwd.h"
32 #include "FWCore/Framework/interface/EDAnalyzer.h"
33
34 #include "FWCore/Framework/interface/Event.h"
35 #include "FWCore/Framework/interface/MakerMacros.h"
36
37 #include "FWCore/ParameterSet/interface/ParameterSet.h"
38 #include "FWCore/Utilities/interface/InputTag.h"
39 #include "DataFormats/TrackReco/interface/Track.h"
40 #include "DataFormats/TrackReco/interface/TrackFwd.h"
41 #include "FWCore/ServiceRegistry/interface/Service.h"
42 #include "CommonTools/UtilAlgos/interface/TFileService.h"
43 #include "UserCode/PVStudy/interface/PVEffAnalyzer.h"
44 #include "DataFormats/Common/interface/RefToBaseVector.h"
45 #include "DataFormats/Common/interface/RefToBase.h"
46 // Vertex
47 #include "DataFormats/VertexReco/interface/Vertex.h"
48 #include "RecoVertex/VertexPrimitives/interface/TransientVertex.h"
49
50 // simulated vertices,..., add <use name=SimDataFormats/Vertex> and <../Track>
51 #include <SimDataFormats/Vertex/interface/SimVertex.h>
52 #include <SimDataFormats/Vertex/interface/SimVertexContainer.h>
53 #include <SimDataFormats/Track/interface/SimTrack.h>
54 #include <SimDataFormats/Track/interface/SimTrackContainer.h>
55 // TrackingParticle Associators
56 #include "SimTracker/TrackAssociation/interface/TrackAssociatorByHits.h"
57 #include "SimTracker/TrackerHitAssociation/interface/TrackerHitAssociator.h"
58 #include "SimTracker/Records/interface/TrackAssociatorRecord.h"
59 // BeamSpot
60 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
61
62 // For the clusters
63 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
64 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
65 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
66 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
67
68 //root
69 #include <TROOT.h>
70 #include <TF1.h>
71 #include <TString.h>
72 #include <TStyle.h>
73 #include <TPaveStats.h>
74 #include <TPad.h>
75
76
77
78 using namespace std;
79 using namespace edm;
80 using namespace reco;
81
82 typedef math::XYZTLorentzVectorF LorentzVector;
83 typedef math::XYZPoint Point;
84
85 PVEffAnalyzer::PVEffAnalyzer(const edm::ParameterSet& iConfig)
86 :theTrackClusterizer(iConfig.getParameter<edm::ParameterSet>("TkClusParameters")),
87 theTrackFilter(iConfig.getParameter<edm::ParameterSet>("TkFilterParameters"))
88 {
89 //=======================================================================
90 // Get configuration for TrackTupleMaker
91 //=======================================================================
92 simG4_ = iConfig.getParameter<edm::InputTag>( "simG4" );
93 trackCollectionTag_ = iConfig.getUntrackedParameter<edm::InputTag>("trackCollection");
94 splitTrackCollection1Tag_ = iConfig.getUntrackedParameter<edm::InputTag>("splitTrackCollection1");
95 splitTrackCollection2Tag_ = iConfig.getUntrackedParameter<edm::InputTag>("splitTrackCollection2");
96 vertexCollectionTag_ = iConfig.getUntrackedParameter<edm::InputTag>("vertexCollection");
97 splitVertexCollection1Tag_ = iConfig.getUntrackedParameter<edm::InputTag>("splitVertexCollection1");
98 splitVertexCollection2Tag_ = iConfig.getUntrackedParameter<edm::InputTag>("splitVertexCollection2");
99 verbose_ = iConfig.getUntrackedParameter<bool>("verbose",false);
100 realData_ = iConfig.getUntrackedParameter<bool>("realData",false);
101 useTP_ = iConfig.getUntrackedParameter<bool>("useTP",false);
102 useAssociator_ = iConfig.getUntrackedParameter<bool>("useAssociator",false);
103 histoFileName_ = iConfig.getUntrackedParameter<std::string> ("histoFileName");
104 outputfilename_ = iConfig.getUntrackedParameter<string>("OutputFileName");
105 nTrkMin_ = iConfig.getUntrackedParameter<int>("nTrkMin");
106 nTrkMax_ = iConfig.getUntrackedParameter<int>("nTrkMax");
107 zsigncut_ = iConfig.getUntrackedParameter<double>("zsigncut");
108 ptcut_ = iConfig.getUntrackedParameter<double>("ptcut");
109 analyze_ = iConfig.getUntrackedParameter<bool>("analyze",false);
110 saventuple_ = iConfig.getUntrackedParameter<bool>("saventuple",false);
111 bsSrc = iConfig.getParameter< edm::InputTag >("beamSpot");
112 reqCluster_ = iConfig.getUntrackedParameter<bool>("reqCluster",false);
113
114 // Specify the data mode vector
115 if(realData_) datamode.push_back(0);
116 else {
117 datamode.push_back(0);
118 datamode.push_back(1);
119 }
120
121 if(saventuple_) {
122 file_ = TFile::Open(outputfilename_.c_str(),"RECREATE");
123 // pvtxtree_ analyzing the pvtxs ootb
124 pvtxtree_ = new TTree("pvtxtree","pvtxtree");
125 if(!realData_) {
126 pvtxtree_->Branch("nsimPV",&nsimPV_,"nsimPV/I");
127 pvtxtree_->Branch("nsimTrkPV",&nsimTrkPV_,"nsimTrkPV[nsimPV]/I");
128 pvtxtree_->Branch("simx",&simx_,"simx[nsimPV]/D");
129 pvtxtree_->Branch("simy",&simy_,"simy[nsimPV]/D");
130 pvtxtree_->Branch("simz",&simz_,"simz[nsimPV]/D");
131 pvtxtree_->Branch("simptsq",&simptsq_,"simptsq[nsimPV]/D");
132 pvtxtree_->Branch("simDeltaZ",&simDeltaZ_,"simDeltaZ[nsimPV]/D");
133 }
134 }
135
136 theFile = new TFile(histoFileName_.c_str(), "RECREATE");
137 theFile->mkdir("Summary");
138 theFile->cd();
139
140
141 // Book MC only plots
142 if (!realData_) {
143 h_gen = new PVEffHistograms();
144 h_gen->Init("generator");
145 }
146
147 h_summary = new PVEffHistograms();
148 //Book histograms sensitive to data/mc
149 for (vector<int>::const_iterator it= datamode.begin(); it != datamode.end() ; ++it) {
150 string suffix;
151 edm::LogInfo("Debug")<<"datamode = "<< *it<<endl;
152 switch(*it) {
153 case 0: suffix = "";
154 break;
155 case 1: suffix = "_mct";
156 break;
157 }
158 h_summary->Init("summary",suffix, nTrkMin_, nTrkMax_);
159 }
160
161 }
162
163 PVEffAnalyzer::~PVEffAnalyzer()
164 {
165 // do anything here that needs to be done at desctruction time
166 // (e.g. close files, deallocate resources etc.)
167 theFile->cd();
168 theFile->cd("Summary");
169 h_summary->Save();
170 if (!realData_)
171 h_gen->Save();
172 theFile->Close();
173 }
174
175 //
176 // member functions
177 //
178 std::vector<PVEffAnalyzer::simPrimaryVertex> PVEffAnalyzer::getSimPVs(const Handle<HepMCProduct> & evtMC, std::string suffix, double & ptcut_)
179 {
180 std::vector<PVEffAnalyzer::simPrimaryVertex> simpv;
181 const HepMC::GenEvent* evt=evtMC->GetEvent();
182 if (evt) {
183 edm::LogInfo("SimPVs") << "[getSimPVs] process id " << evt->signal_process_id()<<endl;
184 edm::LogInfo("SimPVs") << "[getSimPVs] signal process vertex " << ( evt->signal_process_vertex() ?
185 evt->signal_process_vertex()->barcode() : 0 ) <<endl;
186 edm::LogInfo("SimPVs") << "[getSimPVs] number of vertices " << evt->vertices_size() << endl;
187
188 int idx=0; int npv=0;
189 for(HepMC::GenEvent::vertex_const_iterator vitr= evt->vertices_begin();
190 vitr != evt->vertices_end(); ++vitr ) { // loop for vertex ...
191 HepMC::FourVector pos = (*vitr)->position();
192 //HepLorentzVector pos = (*vitr)->position();
193
194 // t component of PV:
195 for ( HepMC::GenVertex::particle_iterator mother = (*vitr)->particles_begin(HepMC::parents);
196 mother != (*vitr)->particles_end(HepMC::parents); ++mother ) {
197 // edm::LogInfo("SimPVs") << "Status = " << (*mother)->status() << endl;
198 HepMC::GenVertex * mv=(*mother)->production_vertex();
199 if( ((*mother)->status() == 3) && (!mv)) {
200 // edm::LogInfo("SimPVs") << "npv= " << npv << endl;
201 if (npv == 0) {
202 h_gen->Fill1d("genPart_cT", pos.t()); // mm
203 h_gen->Fill1d("genPart_T", pos.t()/299.792458); // ns
204 }
205 npv++;
206 }
207 }
208 // if (pos.t()>0) { continue;} // for 22X when t of PV was not smeared
209
210 bool hasMotherVertex=false;
211 if (verbose_) cout << "[getSimPVs] mothers of vertex[" << ++idx << "]: " << endl;
212 for ( HepMC::GenVertex::particle_iterator mother = (*vitr)->particles_begin(HepMC::parents);
213 mother != (*vitr)->particles_end(HepMC::parents); ++mother ) {
214 HepMC::GenVertex * mv=(*mother)->production_vertex();
215 // if (verbose_) cout << "Status = " << (*mother)->status() << endl;
216 if (mv) {
217 hasMotherVertex=true;
218 if(!verbose_) break; //if verbose_, print all particles of gen vertices
219 }
220 if(verbose_) {
221 cout << "\t";
222 (*mother)->print();
223 }
224 }
225
226 if(hasMotherVertex) continue;
227
228 // could be a new vertex, check all primaries found so far to avoid multiple entries
229 const double mm2cm=0.1;
230 simPrimaryVertex sv(pos.x()*mm2cm,pos.y()*mm2cm,pos.z()*mm2cm); // sim unit mm, rec unit cm
231 simPrimaryVertex *vp=NULL; // will become non-NULL if a vertex is found and then point to it
232 for(vector<simPrimaryVertex>::iterator v0=simpv.begin();
233 v0!=simpv.end(); v0++){
234 if( (fabs(sv.x-v0->x)<1e-5) && (fabs(sv.y-v0->y)<1e-5) && (fabs(sv.z-v0->z)<1e-5)){
235 vp=&(*v0);
236 break;
237 }
238 }
239
240 if(!vp){
241 // this is a new vertex
242 edm::LogInfo("SimPVs") << "[getSimPVs] this is a new vertex " << sv.x << " " << sv.y << " " << sv.z << endl;
243 simpv.push_back(sv);
244 vp=&simpv.back();
245 }else{
246 edm::LogInfo("SimPVs") << "[getSimPVs] this is not a new vertex " << sv.x << " " << sv.y << " " << sv.z << endl;
247 }
248 vp->genVertex.push_back((*vitr)->barcode());
249 // collect final state descendants
250 for ( HepMC::GenVertex::particle_iterator daughter = (*vitr)->particles_begin(HepMC::descendants);
251 daughter != (*vitr)->particles_end(HepMC::descendants);
252 ++daughter ) {
253 if (isFinalstateParticle(*daughter)){
254 if ( find(vp->finalstateParticles.begin(), vp->finalstateParticles.end(),(*daughter)->barcode())
255 == vp->finalstateParticles.end()){
256 vp->finalstateParticles.push_back((*daughter)->barcode());
257 HepMC::FourVector m=(*daughter)->momentum();
258 // the next four lines used to be "vp->ptot+=m;" in the days of CLHEP::HepLorentzVector
259 // but adding FourVectors seems not to be foreseen
260 vp->ptot.setPx(vp->ptot.px()+m.px());
261 vp->ptot.setPy(vp->ptot.py()+m.py());
262 vp->ptot.setPz(vp->ptot.pz()+m.pz());
263 vp->ptot.setE(vp->ptot.e()+m.e());
264 vp->ptsq+=(m.perp())*(m.perp());
265 if ( (m.perp()>ptcut_) && (fabs(m.pseudoRapidity())<2.5) && isCharged( *daughter ) ){
266 vp->nGenTrk++;
267 }
268 }
269 }
270 }//loop MC vertices daughters
271 }//loop MC vertices
272 }
273 return simpv;
274 }
275
276
277 // get sim pv from TrackingParticles/TrackingVertex
278 std::vector<PVEffAnalyzer::simPrimaryVertex> PVEffAnalyzer::getSimPVs( const edm::Handle<TrackingVertexCollection> & tVC )
279 {
280 std::vector<simPrimaryVertex> simpv;
281 //std::cout <<"number of vertices " << tVC->size() << std::endl;
282
283 if(verbose_){std::cout << "getSimPVs TrackingVertexCollection " << std::endl;}
284
285 for (TrackingVertexCollection::const_iterator v = tVC -> begin(); v != tVC -> end(); ++v) {
286
287 if(verbose_){
288 std::cout << (v->eventId()).event() << v -> position() << v->g4Vertices().size() <<" " <<v->genVertices().size() << " t=" <<v->position().t()*1.e12 <<" ==0:" <<(v->position().t()>0) << std::endl;
289 for( TrackingVertex::g4v_iterator gv=v->g4Vertices_begin(); gv!=v->g4Vertices_end(); gv++){
290 std::cout << *gv << std::endl;
291 }
292 std::cout << "----------" << std::endl;
293 }
294
295 // bool hasMotherVertex=false;
296 if ((unsigned int) v->eventId().event()<simpv.size()) continue;
297 //if (v->position().t()>0) { continue;} // skip secondary vertices (obsolete + doesn't work)
298 if (fabs(v->position().z())>1000) continue; // skip funny junk vertices
299
300 // could be a new vertex, check all primaries found so far to avoid multiple entries
301 const double mm=1.0; // for tracking vertices
302 simPrimaryVertex sv(v->position().x()*mm,v->position().y()*mm,v->position().z()*mm);
303
304 for (TrackingParticleRefVector::iterator iTrack = v->daughterTracks_begin(); iTrack != v->daughterTracks_end(); ++iTrack){
305 //cout <<((**iTrack).eventId()).event() << " "; // an iterator of Refs, dereference twice. Cool eyh?
306 sv.eventId=(**iTrack).eventId();
307 sv.ptsq += pow((**iTrack).pt(),2);
308 sv.nGenTrk++;
309 }
310
311 simPrimaryVertex *vp=NULL; // will become non-NULL if a vertex is found and then point to it
312 for(std::vector<simPrimaryVertex>::iterator v0=simpv.begin();
313 v0!=simpv.end(); v0++){
314 if( (sv.eventId==v0->eventId) && (fabs(sv.x-v0->x)<1e-5) && (fabs(sv.y-v0->y)<1e-5) && (fabs(sv.z-v0->z)<1e-5)){
315 vp=&(*v0);
316 break;
317 }
318 }
319 if(!vp){
320 // this is a new vertex
321 if(verbose_){std::cout << "this is a new vertex " << sv.eventId.event() << " " << sv.x << " " << sv.y << " " << sv.z <<std::endl;}
322 simpv.push_back(sv);
323 vp=&simpv.back();
324 }else{
325 if(verbose_){std::cout << "this is not a new vertex" << sv.x << " " << sv.y << " " << sv.z <<std::endl;}
326 }
327
328
329 // Loop over daughter track(s)
330 if(verbose_){
331 for (TrackingVertex::tp_iterator iTP = v -> daughterTracks_begin(); iTP != v -> daughterTracks_end(); ++iTP) {
332 std::cout << " Daughter momentum: " << (*(*iTP)).momentum();
333 std::cout << " Daughter type " << (*(*iTP)).pdgId();
334 std::cout << std::endl;
335 }
336 }
337 }
338 if(verbose_){
339 cout << "------- PVEffAnalyzer simPVs from TrackingVertices -------" << endl;
340 for(std::vector<simPrimaryVertex>::iterator v0=simpv.begin();
341 v0!=simpv.end(); v0++){
342 cout << "z=" << v0->z << " event=" << v0->eventId.event() << endl;
343 }
344 cout << "-----------------------------------------------" << endl;
345 }
346 return simpv;
347 }
348
349
350
351
352 // ------------ method called to for each event ------------
353 void
354 PVEffAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
355 {
356 using namespace edm;
357 using namespace reco;
358
359 edm::ESHandle<TransientTrackBuilder> theB;
360 iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
361
362 //========================================================================
363 // Step 0: Prepare root variables and get information from the Event
364 //========================================================================
365
366 edm::LogInfo("Debug")<<"[PVEffAnalyzer]"<<endl;
367
368 // ====== TrackCollection
369 edm::Handle<edm::View<reco::Track> > trackCollectionHandle;
370 iEvent.getByLabel(trackCollectionTag_, trackCollectionHandle);
371
372 // ====== splitTrackCollection1
373 edm::Handle<edm::View<reco::Track> > splitTrackCollection1Handle;
374 iEvent.getByLabel(splitTrackCollection1Tag_, splitTrackCollection1Handle);
375
376 // ====== splitTrackCollection2
377 edm::Handle<edm::View<reco::Track> > splitTrackCollection2Handle;
378 iEvent.getByLabel(splitTrackCollection2Tag_, splitTrackCollection2Handle);
379
380 // ======= PrimaryVertexCollection
381 static const reco::VertexCollection s_empty_vertexColl;
382 const reco::VertexCollection *vertexColl = &s_empty_vertexColl;
383 edm::Handle<reco::VertexCollection> vertexCollectionHandle;
384 iEvent.getByLabel(vertexCollectionTag_, vertexCollectionHandle);
385 if( iEvent.getByLabel(vertexCollectionTag_, vertexCollectionHandle)) {
386 vertexColl = vertexCollectionHandle.product();
387 } else {
388 edm::LogInfo("Debug") << "[PVEffAnalyzer] vertexCollection cannot be found -> using empty collection of same type." <<endl;
389 }
390 // ====== splitVertexCollection1
391 static const reco::VertexCollection s_empty_splitVertexColl1;
392 const reco::VertexCollection *splitVertexColl1 = &s_empty_splitVertexColl1;
393 edm::Handle<reco::VertexCollection> splitVertexCollection1Handle;
394 iEvent.getByLabel(splitVertexCollection1Tag_, splitVertexCollection1Handle);
395 if( iEvent.getByLabel(splitVertexCollection1Tag_, splitVertexCollection1Handle)) {
396 splitVertexColl1 = splitVertexCollection1Handle.product();
397 } else {
398 edm::LogInfo("Debug") << "[PVEffAnalyzer] splitVertexCollection1 cannot be found -> using empty collection of same type." <<endl;
399 }
400 // ====== splitVertexCollection2
401 static const reco::VertexCollection s_empty_splitVertexColl2;
402 const reco::VertexCollection *splitVertexColl2 = &s_empty_splitVertexColl2;
403 edm::Handle<reco::VertexCollection> splitVertexCollection2Handle;
404 iEvent.getByLabel(splitVertexCollection2Tag_, splitVertexCollection2Handle);
405 if( iEvent.getByLabel(splitVertexCollection2Tag_, splitVertexCollection2Handle)) {
406 splitVertexColl2 = splitVertexCollection2Handle.product();
407 } else {
408 edm::LogInfo("Debug") << "[PVEffAnalyzer] splitVertexCollection2 cannot be found -> using empty collection of same type." <<endl;
409 }
410
411
412 // ======== BeamSpot accessors
413 edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
414 iEvent.getByLabel(bsSrc,recoBeamSpotHandle);
415 reco::BeamSpot bs = *recoBeamSpotHandle;
416 const Point beamSpot = recoBeamSpotHandle.isValid() ? Point(recoBeamSpotHandle->x0(), recoBeamSpotHandle->y0(), recoBeamSpotHandle->z0()) : Point(0, 0, 0);
417
418 edm::LogInfo("Debug")<<"[PVEffAnalyzer] End accessing the track, beamSpot, primary vertex collections"<<endl;
419
420 // ========== MC simvtx accessor
421 if (!realData_) {
422 edm::Handle<SimVertexContainer> simVtxs;
423 iEvent.getByLabel( simG4_, simVtxs);
424
425 edm::Handle<SimTrackContainer> simTrks;
426 iEvent.getByLabel( simG4_, simTrks);
427 }
428
429 // ========== GET PDT
430 try{
431 iSetup.getData(pdt);
432 }catch(const Exception&){
433 edm::LogInfo("Debug") << "[PVEffAnalyzer] Some problem occurred with the particle data table. This may not work !." <<endl;
434 }
435
436 //setUpVectors(RealData, nTrkMin_, nTrkMax_ ) ;
437
438 // ======= Analyze MC efficiency
439 if (!realData_) {
440
441 bool MC=false; // the flag to see if the MC level information is complete
442
443 Handle<HepMCProduct> evtMC;
444 iEvent.getByLabel("generator",evtMC);
445
446 edm::Handle<TrackingParticleCollection> TPCollectionH ;
447 bool gotTP=iEvent.getByLabel("mergedtruth","MergedTrackTruth",TPCollectionH);
448
449 edm::Handle<TrackingVertexCollection> TVCollectionH ;
450 bool gotTV=iEvent.getByLabel("mergedtruth","MergedTrackTruth",TVCollectionH);
451
452 std::vector<simPrimaryVertex> simpv;
453
454 if(gotTV){
455 MC=true;
456 if(verbose_){
457 cout << "Found Tracking Vertices " << endl;
458 }
459 simpv=getSimPVs(TVCollectionH);
460 }
461 else if (evtMC.isValid()) {
462 MC=true;
463 edm::LogInfo("Debug") << "[PVEffAnalyzer] generator HepMCProduct found"<< endl;
464 } else {
465 edm::LogInfo("Debug") << "[PVEffAnalyzer] no HepMCProduct found"<< endl;
466 MC=false;
467 }
468
469
470 if(MC){
471 TString suffix = "_mct";
472 // make a list of primary vertices:
473 h_gen->Fill1d("nsimPV", simpv.size());
474
475 int isimpv = 0;
476 nsimPV_ = simpv.size();
477
478 for(int isimpv = 0; isimpv<int(simpv.size());isimpv++) {
479 simPrimaryVertex vsim = simpv.at(isimpv);
480 nsimTrkPV_[isimpv] =vsim.nGenTrk;
481 simx_[isimpv] = vsim.x;
482 simy_[isimpv] = vsim.y;
483 simz_[isimpv] = vsim.z;
484 simptsq_[isimpv] = vsim.ptsq;
485 if(simpv.size() > 1 )
486 if(verbose_)
487 std::cout<<"Simulated Vertex # " << isimpv << ": ptsq = " << vsim.ptsq<<std::endl;
488 if(isimpv >=1 )
489 simDeltaZ_[isimpv-1] = vsim.z-simz_[0];
490 }
491
492 // Just analyze the first vertex in the collection
493 simPrimaryVertex vsim = *simpv.begin();
494 int ntracks(0);
495 if( useTP_ )
496 ntracks = vsim.nGenTrk;
497 else {
498 reco::RecoToSimCollection recSimColl;
499 if (useAssociator_) {
500 // get TrackingParticle
501 edm::Handle<TrackingParticleCollection> TPCollectionHandle ;
502 iEvent.getByLabel("mergedtruth","MergedTrackTruth",TPCollectionHandle);
503 // get associators
504 edm::ESHandle<TrackAssociatorBase> theAssociator;
505 iSetup.get<TrackAssociatorRecord>().get("TrackAssociatorByHitsRecoDenom",theAssociator);
506 recSimColl = (theAssociator.product())->associateRecoToSim( trackCollectionHandle,
507 TPCollectionHandle,
508 & iEvent);
509 }
510
511 // == use the nubmer of tracks that are filtered and clusterized
512 // copy code from CMSSW/RecoVertex/PrimaryVertexProducer/src/PrimaryVertexProducerAlgorithm.cc
513 std::vector<reco::TransientTrack> seltks;
514 for(unsigned int i=0; i<trackCollectionHandle->size(); ++i){
515 edm::RefToBase<reco::Track> track(trackCollectionHandle, i);
516 if(useAssociator_) {
517 bool isSimMatched(false);
518 std::vector<std::pair<TrackingParticleRef, double> > tp;
519 if(recSimColl.find(track) != recSimColl.end()){
520 tp = recSimColl[track];
521 if (tp.size()!=0)
522 isSimMatched = true;
523 }
524 if(!isSimMatched) continue;
525 }
526 TransientTrack t_track = (*theB).build(*track);
527 t_track.setBeamSpot(bs);
528 if (theTrackFilter(t_track)) seltks.push_back(t_track);
529 }
530 std::vector< std::vector<reco::TransientTrack> > clusters = theTrackClusterizer.clusterize(seltks);
531 if (clusters.size() == 1 && (clusters.begin()->size()) >1 ) {
532 for(unsigned idx_tt = 0; idx_tt < clusters.begin()->size(); idx_tt++) {
533 reco::Track track = (*(clusters.begin())).at(idx_tt).track();
534 if(track.pt()>ptcut_) // && TMath::Abs(track.eta())<2.5) // count only good tracks
535 ntracks++;
536 }
537 }
538 } // == end of getting the ntracks
539
540 // == Fill Efficiency Histogram
541 if(ntracks>1) {
542 h_summary->Fill1d(TString("eff_denom_ntrack"+suffix), ntracks);
543 if(!vertexColl->begin()->isFake())
544 h_summary->Fill1d(TString("deltazSign"+suffix), vsim.z - vertexColl->begin()->z());
545 if ( isAssoVertex(vsim, *vertexColl->begin(), zsigncut_) )
546 h_summary->Fill1d(TString("eff_numer_ntrack"+suffix), ntracks);
547 }
548
549 // == Fill FakeRate Histogram
550 if(!vertexColl->begin()->isFake ()) {
551 h_summary->Fill1d(TString("fakerate_denom_ntrack"+suffix), int(vertexColl->begin()->tracksSize()));
552 if (!isAssoVertex(vsim, *vertexColl->begin(), zsigncut_) )
553 h_summary->Fill1d(TString("fakerate_numer_ntrack"+suffix), int(vertexColl->begin()->tracksSize()));
554 }
555 } // end of MC analyzing
556 } // end of MC access
557
558 // == Start of Split-Method
559 // Event selections
560 if ( !isGoodSplitEvent( *vertexColl->begin())) return;
561 for(reco::Vertex::trackRef_iterator t = vertexColl->begin()->tracks_begin();
562 t!=vertexColl->begin()->tracks_end(); t++) {
563 h_summary->Fill1d("trackweight",vertexColl->begin()->trackWeight(*t));
564 }
565
566 // == Start of FakeRate Analyzing
567
568 if(!splitVertexColl1->begin()->isFake()) {
569 h_summary->Fill1d("fakerate_denom_ntrack", int(splitVertexColl1->begin()->tracksSize()));
570 if ( !isGoodProbeVertex( *splitVertexColl1->begin(), *vertexColl->begin(), zsigncut_) )
571 h_summary->Fill1d("fakerate_numer_ntrack", int(splitVertexColl1->begin()->tracksSize()));
572 }
573
574 // == End of FakeRate Analyzing
575 // == Start of Efficiency Analyzing
576 std::vector< std::vector<reco::TransientTrack> > clusters1;
577 std::vector< std::vector<reco::TransientTrack> > clusters2;
578 if(splitTrackCollection1Handle.isValid())
579 clusters1 = theTrackClusterizer.clusterize((*theB).build(splitTrackCollection1Handle));
580 if(splitTrackCollection2Handle.isValid())
581 clusters2 = theTrackClusterizer.clusterize((*theB).build(splitTrackCollection2Handle));
582
583 // clustering requirement
584 if( !goodCluster (clusters1) || !goodCluster (clusters2) ) return;
585
586 // TagVertex Requirement
587 if ( isGoodTagVertex( *(splitVertexColl2->begin()), *(vertexColl->begin()), zsigncut_) ) {
588 if(verbose_) {
589 std::cout<<"splitTrackCollection1Handle->size() = " << int(splitTrackCollection1Handle->size()) << std::endl;
590 std::cout<<"splitTrackCollection2Handle->size() = " << int(splitTrackCollection2Handle->size()) << std::endl;
591 }
592 int nGoodTracks = 0;
593 for(unsigned int i = 0; i < splitTrackCollection1Handle->size(); i++ ) {
594 edm::RefToBase<reco::Track> track(splitTrackCollection1Handle, i);
595 if(track->pt() > ptcut_) // && abs(track->eta()) < 2.5 )
596 nGoodTracks++;
597 }
598
599 h_summary->Fill1d("eff_denom_ntrack", nGoodTracks);
600 // ProbeVertex Requirement
601 if ( isGoodProbeVertex( *splitVertexColl1->begin(), *vertexColl->begin(), zsigncut_) ) {
602 h_summary->Fill1d("eff_numer_ntrack", nGoodTracks);
603 h_summary->Fill1d("avgWeight_orgvtx_eff",avgWeight(*vertexColl->begin()) );
604 h_summary->Fill1d("avgWeight_tagvtx_eff",avgWeight(*splitVertexColl2->begin()) );
605 h_summary->Fill1d("avgWeight_probevtx_eff",avgWeight(*splitVertexColl1->begin()) );
606 }
607 else
608 {
609 if (nGoodTracks >= 2) {
610 h_summary->Fill1d("avgWeight_orgvtx_ineff",avgWeight(*vertexColl->begin()) );
611 h_summary->Fill1d("avgWeight_tagvtx_ineff",avgWeight(*splitVertexColl2->begin()) );
612 h_summary->Fill1d("avgWeight_probevtx_ineff",avgWeight(*splitVertexColl1->begin()) );
613
614 std::cout<<"**********Inefficiency Found************** " <<std::endl;
615 std::cout<<"==== Original Vertex: " <<std::endl;
616 //printRecVtx(*vertexColl->begin());
617 printRecVtxs(vertexCollectionHandle);
618 std::cout<<"==== Probe Vertex: " <<std::endl;
619 //printRecVtx(*splitVertexColl1->begin());
620 printRecVtxs(splitVertexCollection1Handle);
621 std::cout<<"==== Tag Vertex: " <<std::endl;
622 //printRecVtx(*splitVertexColl2->begin());
623 printRecVtxs(splitVertexCollection2Handle);
624 std::cout<<"==== Probe cluster information "<< std::endl;
625 printCluster(clusters1, beamSpot);
626 std::cout<<"==== Tag cluster information "<< std::endl;
627 printCluster(clusters2, beamSpot);
628 std::cout<<"**********Inefficiency Found************** " <<std::endl;
629 }
630 }
631 }
632 // == End of Efficiency Analyzing
633 // == Save Ntuple
634 if(saventuple_) {
635 pvtxtree_->Fill();
636 }
637
638 }
639
640
641 // ------------ method called once each job just before starting event loop ------------
642 void
643 PVEffAnalyzer::beginJob()
644 {
645 }
646
647 // ------------ method called once each job just after ending the event loop ------------
648 void
649 PVEffAnalyzer::endJob() {
650 edm::LogInfo("Analysis") << "[endJob] Analyzing PV info" << endl;
651
652 if(saventuple_) {
653 file_->cd();
654 pvtxtree_->Write();
655 file_->Close();
656 }
657
658
659 for (vector<int>::const_iterator it= datamode.begin(); it != datamode.end() ; ++it) {
660 string suffix;
661 edm::LogInfo("Debug")<<"datamode = "<< *it<<endl;
662 switch(*it) {
663 case 0: suffix = "";
664 break;
665 case 1: suffix = "_mct";
666 break;
667 }
668 if (analyze_) {
669 MakeEff(h_summary->ReadHisto1D(TString("eff_numer_ntrack"+suffix)), h_summary->ReadHisto1D(TString("eff_denom_ntrack"+suffix)), h_summary->ReadHisto1D(TString("eff_ntrack"+suffix)), false, 1);
670 MakeEff(h_summary->ReadHisto1D(TString("fakerate_numer_ntrack"+suffix)), h_summary->ReadHisto1D(TString("fakerate_denom_ntrack"+suffix)), h_summary->ReadHisto1D(TString("fakerate_ntrack"+suffix)), false, 1);
671 }
672 }
673 }
674
675
676
677 bool PVEffAnalyzer::isResonance(const HepMC::GenParticle * p){
678 double ctau=(pdt->particle( abs(p->pdg_id()) ))->lifetime();
679 edm::LogInfo("Debug") << "[isResonance] isResonance " << p->pdg_id() << " " << ctau << endl;
680 return ctau >0 && ctau <1e-6;
681 }
682
683 bool PVEffAnalyzer::isFinalstateParticle(const HepMC::GenParticle * p){
684 return ( !p->end_vertex() && p->status()==1 );
685 }
686
687 bool PVEffAnalyzer::isCharged(const HepMC::GenParticle * p){
688 const ParticleData * part = pdt->particle( p->pdg_id() );
689 if (part){
690 return part->charge()!=0;
691 }else{
692 // the new/improved particle table doesn't know anti-particles
693 return pdt->particle( -p->pdg_id() )!=0;
694 }
695 }
696
697 void PVEffAnalyzer::printSimVtxs(const Handle<SimVertexContainer> simVtxs){
698 int i=0;
699 for(SimVertexContainer::const_iterator vsim=simVtxs->begin();
700 vsim!=simVtxs->end(); ++vsim){
701 cout << i++ << ")" << scientific
702 << " evtid=" << vsim->eventId().event()
703 << " sim x=" << vsim->position().x()
704 << " sim y=" << vsim->position().y()
705 << " sim z=" << vsim->position().z()
706 << " sim t=" << vsim->position().t()
707 << " parent=" << vsim->parentIndex()
708 << endl;
709 }
710 }
711
712 void PVEffAnalyzer::printRecVtxs(const Handle<reco::VertexCollection> recVtxs){
713 int ivtx=0;
714 for(reco::VertexCollection::const_iterator v=recVtxs->begin();
715 v!=recVtxs->end(); ++v){
716 cout << "Recvtx "<< std::setw(3) << std::setfill(' ')<<ivtx++
717 << "#trk " << std::setw(3) << v->tracksSize()
718 << " chi2 " << std::setw(4) << v->chi2()
719 << " ndof " << std::setw(3) << v->ndof() << endl
720 << " x " << std::setw(8) <<std::fixed << std::setprecision(4) << v->x()
721 << " dx " << std::setw(8) << v->xError()<< endl
722 << " y " << std::setw(8) << v->y()
723 << " dy " << std::setw(8) << v->yError()<< endl
724 << " z " << std::setw(8) << v->z()
725 << " dz " << std::setw(8) << v->zError()
726 << endl;
727 }
728 }
729
730 void PVEffAnalyzer::printRecVtx(const reco::Vertex & v){
731
732 cout << "#trk " << std::setw(3) << v.tracksSize()
733 << " chi2 " << std::setw(4) << v.chi2()
734 << " ndof " << std::setw(3) << v.ndof() << endl
735 << " x " << std::setw(8) <<std::fixed << std::setprecision(4) << v.x()
736 << " dx " << std::setw(8) << v.xError()<< endl
737 << " y " << std::setw(8) << v.y()
738 << " dy " << std::setw(8) << v.yError()<< endl
739 << " z " << std::setw(8) << v.z()
740 << " dz " << std::setw(8) << v.zError()
741 << endl;
742 }
743
744 void PVEffAnalyzer::printSimTrks(const Handle<SimTrackContainer> simTrks){
745 cout << " simTrks type, (momentum), vertIndex, genpartIndex" << endl;
746 int i=1;
747 for(SimTrackContainer::const_iterator t=simTrks->begin();
748 t!=simTrks->end(); ++t){
749 //HepMC::GenParticle* gp=evtMC->GetEvent()->particle( (*t).genpartIndex() );
750 cout << i++ << ")"
751 << (*t)
752 << " index="
753 << (*t).genpartIndex();
754 //if (gp) {
755 // HepMC::GenVertex *gv=gp->production_vertex();
756 // cout << " genvertex =" << (*gv);
757 //}
758 cout << endl;
759 }
760 }
761
762 // Select based on the highest SumPtSq Vertex
763 // Require at least 6 tracks with track weight > ~ 0.8
764 // vtx.ndof = 2*Sum(trackWeight) - 3
765
766 bool PVEffAnalyzer::isGoodSplitEvent( const reco::Vertex & org_vtx)
767 {
768 if ( org_vtx.tracksSize() >= 6 && org_vtx.ndof() >= 6 && avgWeight(org_vtx) > 0.75 ) return true;
769 //if ( org_vtx.tracksSize() >= 6 && org_vtx.ndof() >= 6 ) return true;
770 else
771 return false;
772 }
773
774 // Comparing the tag vertex with the unsplit vertex position in z
775 bool PVEffAnalyzer::isGoodTagVertex( const reco::Vertex & tag_vtx, const reco::Vertex & org_vtx, double & zsign_cut)
776 {
777 if(tag_vtx.isValid () && !tag_vtx.isFake()) {
778 float zsign = TMath::Abs(tag_vtx.z() - org_vtx.z() ) / TMath::Max(tag_vtx.zError(), org_vtx.zError() );
779 //float zsign = TMath::Abs(tag_vtx.z() - org_vtx.z() ) / sqrt(pow(tag_vtx.zError(),2)+pow(org_vtx.zError(),2) );
780 return (zsign < zsign_cut) ;
781 }
782 else
783 return false;
784 }
785
786 // Comparing the probe vertex with the unsplit vertex position in z
787 bool PVEffAnalyzer::isGoodProbeVertex( const reco::Vertex & probe_vtx, const reco::Vertex & org_vtx, double & zsign_cut)
788 {
789 if(probe_vtx.isValid () && !probe_vtx.isFake()) {
790 float zsign = TMath::Abs(probe_vtx.z() - org_vtx.z() ) / TMath::Max(probe_vtx.zError(), org_vtx.zError() );
791 //float zsign = TMath::Abs(probe_vtx.z() - org_vtx.z() ) / sqrt(pow(probe_vtx.zError(),2)+pow(org_vtx.zError(),2) );
792 return (zsign < zsign_cut) ;
793 }
794 else
795 return false;
796 }
797
798
799 bool PVEffAnalyzer::isAssoVertex(const PVEffAnalyzer::simPrimaryVertex & vsim,
800 const reco::Vertex & vrec, double & zsign_cut)
801 {
802 if(vrec.isValid() && !vrec.isFake() )
803 return (TMath::Abs(vsim.z-vrec.z()) < zsign_cut * vrec.zError() );
804 else
805 return false;
806 }
807
808
809
810 void PVEffAnalyzer::MakeEff(TH1D* numer, TH1D* denom, TH1D* eff, //TGraphAsymmErrors* & gr_eff,
811 const bool rebin, const Float_t n ) {
812 eff->Divide( numer, denom, 1,1,"B" );
813
814 if( rebin ) {
815 std::vector<Double_t> binEdges;
816 Int_t bin = denom->GetNbinsX() + 1;
817 binEdges.push_back(denom->GetBinLowEdge(bin));
818 Float_t nEntries = 0;
819 while (bin > 1) {
820 bin --;
821 nEntries = denom->GetBinContent(bin);
822 while (nEntries < n && bin > 1) {
823 bin --;
824 nEntries += denom->GetBinContent(bin);
825 }
826 binEdges.push_back(denom->GetBinLowEdge(bin));
827 }
828
829 Double_t *array = new Double_t[binEdges.size()];
830
831 Int_t j = 0;
832 for (Int_t i = binEdges.size(); i > 0; --i) {
833 array[j] = binEdges[i - 1];
834 ++j;
835 }
836 }
837
838 for (int i = 1; i<eff->GetNbinsX()+1; i++) {
839 if(numer->GetBinContent(i) == 0 || denom->GetBinContent(i) == 0 ) continue;
840 float error = sqrt(eff->GetBinContent(i)*(1-eff->GetBinContent(i))/denom->GetBinContent(i)); // Binominal-Error
841 eff->SetBinError(i, error);
842 }
843 }
844
845 bool PVEffAnalyzer::goodCluster(std::vector< std::vector<reco::TransientTrack> > & clusters)
846 {
847 if( clusters.size() == 1 && ( clusters.begin()->size() > 1 )) return true;
848 else
849 return false;
850 }
851
852 void PVEffAnalyzer::printCluster(std::vector< std::vector<reco::TransientTrack> > & clusters, const Point & beamSpot)
853 {
854 //std::cout<<"theTrackClusterizer.zSeparation()"<<theTrackClusterizer.zSeparation()<<std::endl;
855 cout << "#clusters " << std::setw(3) << clusters.size() << endl;
856 int idx_cluster = 0;
857 for (vector< vector<reco::TransientTrack> >::const_iterator iclus = clusters.begin(); iclus != clusters.end(); iclus++, idx_cluster++) {
858 if((*iclus).size() == 0) continue;
859 cout << "=== Cluster # " << std::setw(3) << idx_cluster << std::setw(3) << " nTracks" << std::setw(3) << (*iclus).size() <<endl;
860 for (int idx_track = 0; idx_track < int((*iclus).size()) ; idx_track++) {
861 reco::Track t = (*iclus).at(idx_track).track();
862 cout << "== Trk # " << std::setw(4) << idx_track
863 << " nHit " << std::setw(4) << t.hitPattern().numberOfValidHits() <<endl
864 << " pt " << std::setw(9) << std::fixed << std::setprecision(4) << t.pt()
865 << " dpt " << std::setw(9) << t.ptError() << endl
866 << " dxy(BS) " << std::setw(9) << t.dxy(beamSpot)
867 << " sigmaDxy " << std::setw(9) << t.dxyError() << endl
868 << " dz(BS) " << std::setw(9) << t.dz(beamSpot)
869 << " sigmaDz " << std::setw(9) << t.dzError() << endl
870 << " xPCA " << std::setw(9) << t.vertex().x()
871 << " yPCA " << std::setw(9) << t.vertex().y()
872 << " zPCA " << std::setw(9) << t.vertex().z()
873 << endl;
874 }
875 }
876
877 }
878 double PVEffAnalyzer::avgWeight(const reco::Vertex & vtx)
879 {
880 if(!vtx.isValid() || vtx.isFake() )
881 return 0;
882 else
883 return (vtx.ndof()+3.0)/2.0/double(vtx.tracksSize());
884 }
885
886
887 void PVEffAnalyzer::SetVarToZero() {
888 nsimPV_ = 0;
889 //simpv
890 for(int i = 0; i < nMaxPVs_; i++) {
891 nsimTrkPV_[i] = 0;
892 simx_[i] = 0;
893 simy_[i] = 0;
894 simz_[i] = 0;
895 simptsq_[i] = 0;
896 if(i>=1)
897 simDeltaZ_[i-1] = 0;
898 }
899 }