1 |
yygao |
1.2 |
|
2 |
yygao |
1.1 |
// -*- 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 |
yygao |
1.5 |
// $Id: PVEffAnalyzer.cc,v 1.4 2010/04/20 12:39:48 yygao Exp $
|
19 |
yygao |
1.1 |
//
|
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 |
yygao |
1.5 |
#include "DataFormats/Common/interface/RefToBaseVector.h"
|
45 |
|
|
#include "DataFormats/Common/interface/RefToBase.h"
|
46 |
|
|
// Vertex
|
47 |
yygao |
1.1 |
#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 |
yygao |
1.5 |
// TrackingParticle Associators
|
56 |
|
|
#include "SimTracker/TrackAssociation/interface/TrackAssociatorByHits.h"
|
57 |
|
|
#include "SimTracker/TrackerHitAssociation/interface/TrackerHitAssociator.h"
|
58 |
|
|
#include "SimTracker/Records/interface/TrackAssociatorRecord.h"
|
59 |
yygao |
1.1 |
// BeamSpot
|
60 |
|
|
#include "DataFormats/BeamSpot/interface/BeamSpot.h"
|
61 |
yygao |
1.5 |
|
62 |
yygao |
1.2 |
// 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 |
yygao |
1.1 |
|
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 |
yygao |
1.5 |
|
77 |
|
|
|
78 |
yygao |
1.1 |
using namespace std;
|
79 |
yygao |
1.2 |
using namespace edm;
|
80 |
|
|
using namespace reco;
|
81 |
|
|
|
82 |
yygao |
1.1 |
typedef math::XYZTLorentzVectorF LorentzVector;
|
83 |
|
|
typedef math::XYZPoint Point;
|
84 |
|
|
|
85 |
|
|
PVEffAnalyzer::PVEffAnalyzer(const edm::ParameterSet& iConfig)
|
86 |
yygao |
1.5 |
:theTrackClusterizer(iConfig.getParameter<edm::ParameterSet>("TkClusParameters")),
|
87 |
|
|
theTrackFilter(iConfig.getParameter<edm::ParameterSet>("TkFilterParameters"))
|
88 |
yygao |
1.1 |
{
|
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 |
yygao |
1.5 |
useTP_ = iConfig.getUntrackedParameter<bool>("useTP",false);
|
102 |
|
|
useAssociator_ = iConfig.getUntrackedParameter<bool>("useAssociator",false);
|
103 |
yygao |
1.1 |
histoFileName_ = iConfig.getUntrackedParameter<std::string> ("histoFileName");
|
104 |
|
|
nTrkMin_ = iConfig.getUntrackedParameter<int>("nTrkMin");
|
105 |
|
|
nTrkMax_ = iConfig.getUntrackedParameter<int>("nTrkMax");
|
106 |
|
|
zsigncut_ = iConfig.getUntrackedParameter<double>("zsigncut");
|
107 |
yygao |
1.4 |
ptcut_ = iConfig.getUntrackedParameter<double>("ptcut");
|
108 |
yygao |
1.2 |
analyze_ = iConfig.getUntrackedParameter<bool>("analyze",false);
|
109 |
yygao |
1.1 |
bsSrc = iConfig.getParameter< edm::InputTag >("beamSpot");
|
110 |
yygao |
1.2 |
reqCluster_ = iConfig.getUntrackedParameter<bool>("reqCluster",false);
|
111 |
|
|
|
112 |
yygao |
1.1 |
// Specify the data mode vector
|
113 |
|
|
if(realData_) datamode.push_back(0);
|
114 |
|
|
else {
|
115 |
|
|
datamode.push_back(0);
|
116 |
|
|
datamode.push_back(1);
|
117 |
|
|
}
|
118 |
|
|
|
119 |
|
|
theFile = new TFile(histoFileName_.c_str(), "RECREATE");
|
120 |
|
|
theFile->mkdir("Summary");
|
121 |
|
|
theFile->cd();
|
122 |
|
|
|
123 |
|
|
// Book MC only plots
|
124 |
|
|
if (!realData_) {
|
125 |
|
|
h_gen = new PVEffHistograms();
|
126 |
|
|
h_gen->Init("generator");
|
127 |
|
|
}
|
128 |
|
|
|
129 |
|
|
h_summary = new PVEffHistograms();
|
130 |
|
|
//Book histograms sensitive to data/mc
|
131 |
|
|
for (vector<int>::const_iterator it= datamode.begin(); it != datamode.end() ; ++it) {
|
132 |
|
|
string suffix;
|
133 |
|
|
edm::LogInfo("Debug")<<"datamode = "<< *it<<endl;
|
134 |
|
|
switch(*it) {
|
135 |
|
|
case 0: suffix = "";
|
136 |
|
|
break;
|
137 |
|
|
case 1: suffix = "_mct";
|
138 |
|
|
break;
|
139 |
|
|
}
|
140 |
|
|
h_summary->Init("summary",suffix, nTrkMin_, nTrkMax_);
|
141 |
|
|
}
|
142 |
|
|
|
143 |
|
|
}
|
144 |
|
|
|
145 |
|
|
PVEffAnalyzer::~PVEffAnalyzer()
|
146 |
|
|
{
|
147 |
|
|
// do anything here that needs to be done at desctruction time
|
148 |
|
|
// (e.g. close files, deallocate resources etc.)
|
149 |
|
|
theFile->cd();
|
150 |
|
|
theFile->cd("Summary");
|
151 |
|
|
h_summary->Save();
|
152 |
|
|
if (!realData_)
|
153 |
|
|
h_gen->Save();
|
154 |
|
|
theFile->Close();
|
155 |
|
|
}
|
156 |
|
|
|
157 |
|
|
//
|
158 |
|
|
// member functions
|
159 |
|
|
//
|
160 |
yygao |
1.4 |
std::vector<PVEffAnalyzer::simPrimaryVertex> PVEffAnalyzer::getSimPVs(const Handle<HepMCProduct> & evtMC, std::string suffix, double & ptcut_)
|
161 |
yygao |
1.1 |
{
|
162 |
|
|
std::vector<PVEffAnalyzer::simPrimaryVertex> simpv;
|
163 |
|
|
const HepMC::GenEvent* evt=evtMC->GetEvent();
|
164 |
|
|
if (evt) {
|
165 |
|
|
edm::LogInfo("SimPVs") << "[getSimPVs] process id " << evt->signal_process_id()<<endl;
|
166 |
|
|
edm::LogInfo("SimPVs") << "[getSimPVs] signal process vertex " << ( evt->signal_process_vertex() ?
|
167 |
|
|
evt->signal_process_vertex()->barcode() : 0 ) <<endl;
|
168 |
|
|
edm::LogInfo("SimPVs") << "[getSimPVs] number of vertices " << evt->vertices_size() << endl;
|
169 |
|
|
|
170 |
|
|
int idx=0; int npv=0;
|
171 |
|
|
for(HepMC::GenEvent::vertex_const_iterator vitr= evt->vertices_begin();
|
172 |
|
|
vitr != evt->vertices_end(); ++vitr ) { // loop for vertex ...
|
173 |
|
|
HepMC::FourVector pos = (*vitr)->position();
|
174 |
|
|
//HepLorentzVector pos = (*vitr)->position();
|
175 |
|
|
|
176 |
|
|
// t component of PV:
|
177 |
|
|
for ( HepMC::GenVertex::particle_iterator mother = (*vitr)->particles_begin(HepMC::parents);
|
178 |
|
|
mother != (*vitr)->particles_end(HepMC::parents); ++mother ) {
|
179 |
|
|
// edm::LogInfo("SimPVs") << "Status = " << (*mother)->status() << endl;
|
180 |
|
|
HepMC::GenVertex * mv=(*mother)->production_vertex();
|
181 |
|
|
if( ((*mother)->status() == 3) && (!mv)) {
|
182 |
|
|
// edm::LogInfo("SimPVs") << "npv= " << npv << endl;
|
183 |
|
|
if (npv == 0) {
|
184 |
|
|
h_gen->Fill1d("genPart_cT", pos.t()); // mm
|
185 |
|
|
h_gen->Fill1d("genPart_T", pos.t()/299.792458); // ns
|
186 |
|
|
}
|
187 |
|
|
npv++;
|
188 |
|
|
}
|
189 |
|
|
}
|
190 |
|
|
// if (pos.t()>0) { continue;} // for 22X when t of PV was not smeared
|
191 |
|
|
|
192 |
|
|
bool hasMotherVertex=false;
|
193 |
|
|
if (verbose_) cout << "[getSimPVs] mothers of vertex[" << ++idx << "]: " << endl;
|
194 |
|
|
for ( HepMC::GenVertex::particle_iterator mother = (*vitr)->particles_begin(HepMC::parents);
|
195 |
|
|
mother != (*vitr)->particles_end(HepMC::parents); ++mother ) {
|
196 |
|
|
HepMC::GenVertex * mv=(*mother)->production_vertex();
|
197 |
|
|
// if (verbose_) cout << "Status = " << (*mother)->status() << endl;
|
198 |
|
|
if (mv) {
|
199 |
|
|
hasMotherVertex=true;
|
200 |
|
|
if(!verbose_) break; //if verbose_, print all particles of gen vertices
|
201 |
|
|
}
|
202 |
|
|
if(verbose_) {
|
203 |
|
|
cout << "\t";
|
204 |
|
|
(*mother)->print();
|
205 |
|
|
}
|
206 |
|
|
}
|
207 |
|
|
|
208 |
|
|
if(hasMotherVertex) continue;
|
209 |
|
|
|
210 |
|
|
// could be a new vertex, check all primaries found so far to avoid multiple entries
|
211 |
|
|
const double mm2cm=0.1;
|
212 |
|
|
simPrimaryVertex sv(pos.x()*mm2cm,pos.y()*mm2cm,pos.z()*mm2cm); // sim unit mm, rec unit cm
|
213 |
|
|
simPrimaryVertex *vp=NULL; // will become non-NULL if a vertex is found and then point to it
|
214 |
|
|
for(vector<simPrimaryVertex>::iterator v0=simpv.begin();
|
215 |
|
|
v0!=simpv.end(); v0++){
|
216 |
|
|
if( (fabs(sv.x-v0->x)<1e-5) && (fabs(sv.y-v0->y)<1e-5) && (fabs(sv.z-v0->z)<1e-5)){
|
217 |
|
|
vp=&(*v0);
|
218 |
|
|
break;
|
219 |
|
|
}
|
220 |
|
|
}
|
221 |
|
|
|
222 |
|
|
if(!vp){
|
223 |
|
|
// this is a new vertex
|
224 |
|
|
edm::LogInfo("SimPVs") << "[getSimPVs] this is a new vertex " << sv.x << " " << sv.y << " " << sv.z << endl;
|
225 |
|
|
simpv.push_back(sv);
|
226 |
|
|
vp=&simpv.back();
|
227 |
|
|
}else{
|
228 |
|
|
edm::LogInfo("SimPVs") << "[getSimPVs] this is not a new vertex " << sv.x << " " << sv.y << " " << sv.z << endl;
|
229 |
|
|
}
|
230 |
|
|
vp->genVertex.push_back((*vitr)->barcode());
|
231 |
|
|
// collect final state descendants
|
232 |
|
|
for ( HepMC::GenVertex::particle_iterator daughter = (*vitr)->particles_begin(HepMC::descendants);
|
233 |
|
|
daughter != (*vitr)->particles_end(HepMC::descendants);
|
234 |
|
|
++daughter ) {
|
235 |
|
|
if (isFinalstateParticle(*daughter)){
|
236 |
|
|
if ( find(vp->finalstateParticles.begin(), vp->finalstateParticles.end(),(*daughter)->barcode())
|
237 |
|
|
== vp->finalstateParticles.end()){
|
238 |
|
|
vp->finalstateParticles.push_back((*daughter)->barcode());
|
239 |
|
|
HepMC::FourVector m=(*daughter)->momentum();
|
240 |
|
|
// the next four lines used to be "vp->ptot+=m;" in the days of CLHEP::HepLorentzVector
|
241 |
|
|
// but adding FourVectors seems not to be foreseen
|
242 |
|
|
vp->ptot.setPx(vp->ptot.px()+m.px());
|
243 |
|
|
vp->ptot.setPy(vp->ptot.py()+m.py());
|
244 |
|
|
vp->ptot.setPz(vp->ptot.pz()+m.pz());
|
245 |
|
|
vp->ptot.setE(vp->ptot.e()+m.e());
|
246 |
|
|
vp->ptsq+=(m.perp())*(m.perp());
|
247 |
yygao |
1.4 |
if ( (m.perp()>ptcut_) && (fabs(m.pseudoRapidity())<2.5) && isCharged( *daughter ) ){
|
248 |
yygao |
1.1 |
vp->nGenTrk++;
|
249 |
|
|
}
|
250 |
|
|
}
|
251 |
|
|
}
|
252 |
|
|
}//loop MC vertices daughters
|
253 |
|
|
}//loop MC vertices
|
254 |
|
|
}
|
255 |
|
|
return simpv;
|
256 |
|
|
}
|
257 |
|
|
|
258 |
|
|
// ------------ method called to for each event ------------
|
259 |
|
|
void
|
260 |
|
|
PVEffAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
|
261 |
|
|
{
|
262 |
|
|
using namespace edm;
|
263 |
|
|
using namespace reco;
|
264 |
yygao |
1.5 |
|
265 |
|
|
edm::ESHandle<TransientTrackBuilder> theB;
|
266 |
|
|
iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
|
267 |
|
|
|
268 |
yygao |
1.1 |
//========================================================================
|
269 |
|
|
// Step 0: Prepare root variables and get information from the Event
|
270 |
|
|
//========================================================================
|
271 |
|
|
|
272 |
|
|
edm::LogInfo("Debug")<<"[PVEffAnalyzer]"<<endl;
|
273 |
|
|
|
274 |
|
|
// ====== TrackCollection
|
275 |
yygao |
1.5 |
edm::Handle<edm::View<reco::Track> > trackCollectionHandle;
|
276 |
yygao |
1.1 |
iEvent.getByLabel(trackCollectionTag_, trackCollectionHandle);
|
277 |
yygao |
1.5 |
|
278 |
yygao |
1.1 |
// ====== splitTrackCollection1
|
279 |
yygao |
1.5 |
edm::Handle<edm::View<reco::Track> > splitTrackCollection1Handle;
|
280 |
yygao |
1.1 |
iEvent.getByLabel(splitTrackCollection1Tag_, splitTrackCollection1Handle);
|
281 |
yygao |
1.5 |
|
282 |
yygao |
1.1 |
// ====== splitTrackCollection2
|
283 |
yygao |
1.5 |
edm::Handle<edm::View<reco::Track> > splitTrackCollection2Handle;
|
284 |
yygao |
1.1 |
iEvent.getByLabel(splitTrackCollection2Tag_, splitTrackCollection2Handle);
|
285 |
|
|
|
286 |
|
|
// ======= PrimaryVertexCollection
|
287 |
|
|
static const reco::VertexCollection s_empty_vertexColl;
|
288 |
|
|
const reco::VertexCollection *vertexColl = &s_empty_vertexColl;
|
289 |
|
|
edm::Handle<reco::VertexCollection> vertexCollectionHandle;
|
290 |
|
|
iEvent.getByLabel(vertexCollectionTag_, vertexCollectionHandle);
|
291 |
|
|
if( iEvent.getByLabel(vertexCollectionTag_, vertexCollectionHandle)) {
|
292 |
|
|
vertexColl = vertexCollectionHandle.product();
|
293 |
|
|
} else {
|
294 |
|
|
edm::LogInfo("Debug") << "[PVEffAnalyzer] vertexCollection cannot be found -> using empty collection of same type." <<endl;
|
295 |
|
|
}
|
296 |
|
|
// ====== splitVertexCollection1
|
297 |
|
|
static const reco::VertexCollection s_empty_splitVertexColl1;
|
298 |
|
|
const reco::VertexCollection *splitVertexColl1 = &s_empty_splitVertexColl1;
|
299 |
|
|
edm::Handle<reco::VertexCollection> splitVertexCollection1Handle;
|
300 |
|
|
iEvent.getByLabel(splitVertexCollection1Tag_, splitVertexCollection1Handle);
|
301 |
|
|
if( iEvent.getByLabel(splitVertexCollection1Tag_, splitVertexCollection1Handle)) {
|
302 |
|
|
splitVertexColl1 = splitVertexCollection1Handle.product();
|
303 |
|
|
} else {
|
304 |
|
|
edm::LogInfo("Debug") << "[PVEffAnalyzer] splitVertexCollection1 cannot be found -> using empty collection of same type." <<endl;
|
305 |
|
|
}
|
306 |
|
|
// ====== splitVertexCollection2
|
307 |
|
|
static const reco::VertexCollection s_empty_splitVertexColl2;
|
308 |
|
|
const reco::VertexCollection *splitVertexColl2 = &s_empty_splitVertexColl2;
|
309 |
|
|
edm::Handle<reco::VertexCollection> splitVertexCollection2Handle;
|
310 |
|
|
iEvent.getByLabel(splitVertexCollection2Tag_, splitVertexCollection2Handle);
|
311 |
|
|
if( iEvent.getByLabel(splitVertexCollection2Tag_, splitVertexCollection2Handle)) {
|
312 |
|
|
splitVertexColl2 = splitVertexCollection2Handle.product();
|
313 |
|
|
} else {
|
314 |
|
|
edm::LogInfo("Debug") << "[PVEffAnalyzer] splitVertexCollection2 cannot be found -> using empty collection of same type." <<endl;
|
315 |
|
|
}
|
316 |
|
|
|
317 |
|
|
|
318 |
|
|
// ======== BeamSpot accessors
|
319 |
|
|
edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
|
320 |
|
|
iEvent.getByLabel(bsSrc,recoBeamSpotHandle);
|
321 |
|
|
reco::BeamSpot bs = *recoBeamSpotHandle;
|
322 |
|
|
const Point beamSpot = recoBeamSpotHandle.isValid() ? Point(recoBeamSpotHandle->x0(), recoBeamSpotHandle->y0(), recoBeamSpotHandle->z0()) : Point(0, 0, 0);
|
323 |
|
|
|
324 |
|
|
edm::LogInfo("Debug")<<"[PVEffAnalyzer] End accessing the track, beamSpot, primary vertex collections"<<endl;
|
325 |
|
|
|
326 |
|
|
// ========== MC simvtx accessor
|
327 |
|
|
if (!realData_) {
|
328 |
|
|
edm::Handle<SimVertexContainer> simVtxs;
|
329 |
|
|
iEvent.getByLabel( simG4_, simVtxs);
|
330 |
|
|
|
331 |
|
|
edm::Handle<SimTrackContainer> simTrks;
|
332 |
|
|
iEvent.getByLabel( simG4_, simTrks);
|
333 |
|
|
}
|
334 |
|
|
|
335 |
|
|
// ========== GET PDT
|
336 |
|
|
try{
|
337 |
|
|
iSetup.getData(pdt);
|
338 |
|
|
}catch(const Exception&){
|
339 |
|
|
edm::LogInfo("Debug") << "[PVEffAnalyzer] Some problem occurred with the particle data table. This may not work !." <<endl;
|
340 |
|
|
}
|
341 |
|
|
|
342 |
|
|
//setUpVectors(RealData, nTrkMin_, nTrkMax_ ) ;
|
343 |
|
|
|
344 |
|
|
// ======= Analyze MC efficiency
|
345 |
|
|
if (!realData_) {
|
346 |
|
|
bool MC=false;
|
347 |
|
|
Handle<HepMCProduct> evtMC;
|
348 |
|
|
iEvent.getByLabel("generator",evtMC);
|
349 |
|
|
if (!evtMC.isValid()) {
|
350 |
|
|
MC=false;
|
351 |
|
|
edm::LogInfo("Debug") << "[PVEffAnalyzer] no HepMCProduct found"<< endl;
|
352 |
|
|
} else {
|
353 |
|
|
edm::LogInfo("Debug") << "[PVEffAnalyzer] generator HepMCProduct found"<< endl;
|
354 |
|
|
MC=true;
|
355 |
|
|
}
|
356 |
|
|
if(MC){
|
357 |
|
|
TString suffix = "_mct";
|
358 |
|
|
// make a list of primary vertices:
|
359 |
|
|
std::vector<simPrimaryVertex> simpv;
|
360 |
yygao |
1.4 |
simpv=getSimPVs(evtMC,"",ptcut_);
|
361 |
yygao |
1.1 |
h_gen->Fill1d("nsimPV", simpv.size());
|
362 |
|
|
|
363 |
|
|
int isimpv = 0;
|
364 |
|
|
|
365 |
|
|
for(std::vector<simPrimaryVertex>::iterator vsim=simpv.begin();
|
366 |
|
|
vsim!=simpv.end(); vsim++, isimpv++){
|
367 |
|
|
//nsimTrkPV_[isimpv] =vsim->nGenTrk;
|
368 |
|
|
//simx_[isimpv] = vsim->x;
|
369 |
|
|
//simy_[isimpv] = vsim->y;
|
370 |
|
|
//simz_[isimpv] = vsim->y;
|
371 |
|
|
//simptsq_[isimpv] = vsim->ptsq;
|
372 |
|
|
if(verbose_ && simpv.size() > 1 )
|
373 |
|
|
std::cout<<"Simulated Vertex # " << isimpv << ": ptsq = " << vsim->ptsq<<std::endl;
|
374 |
|
|
}
|
375 |
|
|
|
376 |
|
|
// Just analyze the first vertex in the collection
|
377 |
|
|
simPrimaryVertex vsim = *simpv.begin();
|
378 |
yygao |
1.5 |
int ntracks(0);
|
379 |
|
|
if( useTP_ )
|
380 |
|
|
ntracks = vsim.nGenTrk;
|
381 |
|
|
else {
|
382 |
|
|
reco::RecoToSimCollection recSimColl;
|
383 |
|
|
if (useAssociator_) {
|
384 |
|
|
// get TrackingParticle
|
385 |
|
|
edm::Handle<TrackingParticleCollection> TPCollectionHandle ;
|
386 |
|
|
iEvent.getByLabel("mergedtruth","MergedTrackTruth",TPCollectionHandle);
|
387 |
|
|
// get associators
|
388 |
|
|
edm::ESHandle<TrackAssociatorBase> theAssociator;
|
389 |
|
|
iSetup.get<TrackAssociatorRecord>().get("TrackAssociatorByHitsRecoDenom",theAssociator);
|
390 |
|
|
recSimColl = (theAssociator.product())->associateRecoToSim( trackCollectionHandle,
|
391 |
|
|
TPCollectionHandle,
|
392 |
|
|
& iEvent);
|
393 |
|
|
}
|
394 |
|
|
|
395 |
|
|
// == use the nubmer of tracks that are filtered and clusterized
|
396 |
|
|
// copy code from CMSSW/RecoVertex/PrimaryVertexProducer/src/PrimaryVertexProducerAlgorithm.cc
|
397 |
|
|
std::vector<reco::TransientTrack> seltks;
|
398 |
|
|
for(unsigned int i=0; i<trackCollectionHandle->size(); ++i){
|
399 |
|
|
edm::RefToBase<reco::Track> track(trackCollectionHandle, i);
|
400 |
|
|
if(useAssociator_) {
|
401 |
|
|
bool isSimMatched(false);
|
402 |
|
|
std::vector<std::pair<TrackingParticleRef, double> > tp;
|
403 |
|
|
if(recSimColl.find(track) != recSimColl.end()){
|
404 |
|
|
tp = recSimColl[track];
|
405 |
|
|
if (tp.size()!=0)
|
406 |
|
|
isSimMatched = true;
|
407 |
|
|
}
|
408 |
|
|
if(!isSimMatched) continue;
|
409 |
|
|
}
|
410 |
|
|
TransientTrack t_track = (*theB).build(*track);
|
411 |
|
|
t_track.setBeamSpot(bs);
|
412 |
|
|
if (theTrackFilter(t_track)) seltks.push_back(t_track);
|
413 |
|
|
}
|
414 |
|
|
std::vector< std::vector<reco::TransientTrack> > clusters = theTrackClusterizer.clusterize(seltks);
|
415 |
|
|
if (clusters.size() == 1 && (clusters.begin()->size()) >1 ) {
|
416 |
|
|
for(unsigned idx_tt = 0; idx_tt < clusters.begin()->size(); idx_tt++) {
|
417 |
|
|
reco::Track track = (*(clusters.begin())).at(idx_tt).track();
|
418 |
|
|
if(track.pt()>ptcut_) // && TMath::Abs(track.eta())<2.5) // count only good tracks
|
419 |
|
|
ntracks++;
|
420 |
|
|
}
|
421 |
|
|
}
|
422 |
|
|
}
|
423 |
|
|
|
424 |
|
|
// == end of getting the ntracks
|
425 |
|
|
if(ntracks>1) {
|
426 |
|
|
h_summary->Fill1d(TString("denom_ntrack"+suffix), ntracks);
|
427 |
|
|
if(!vertexColl->begin()->isFake())
|
428 |
|
|
h_summary->Fill1d(TString("deltazSign"+suffix), vsim.z - vertexColl->begin()->z());
|
429 |
|
|
if ( isAssoVertex(vsim, *vertexColl->begin(), zsigncut_) )
|
430 |
|
|
h_summary->Fill1d(TString("numer_ntrack"+suffix), ntracks);
|
431 |
|
|
}
|
432 |
yygao |
1.1 |
}
|
433 |
|
|
} // End of Analyzing MC Efficiency
|
434 |
|
|
|
435 |
yygao |
1.5 |
// == Start of Split-Method
|
436 |
|
|
// Event selections
|
437 |
|
|
if ( !isGoodSplitEvent( *vertexColl->begin())) return;
|
438 |
|
|
for(reco::Vertex::trackRef_iterator t = vertexColl->begin()->tracks_begin();
|
439 |
|
|
t!=vertexColl->begin()->tracks_end(); t++) {
|
440 |
|
|
h_summary->Fill1d("trackweight",vertexColl->begin()->trackWeight(*t));
|
441 |
|
|
}
|
442 |
|
|
std::vector< std::vector<reco::TransientTrack> > clusters1;
|
443 |
|
|
std::vector< std::vector<reco::TransientTrack> > clusters2;
|
444 |
|
|
if(splitTrackCollection1Handle.isValid())
|
445 |
|
|
clusters1 = theTrackClusterizer.clusterize((*theB).build(splitTrackCollection1Handle));
|
446 |
|
|
if(splitTrackCollection2Handle.isValid())
|
447 |
|
|
clusters2 = theTrackClusterizer.clusterize((*theB).build(splitTrackCollection2Handle));
|
448 |
yygao |
1.2 |
|
449 |
yygao |
1.5 |
// clustering requirement
|
450 |
|
|
if( !goodCluster (clusters1) || !goodCluster (clusters2) ) return;
|
451 |
yygao |
1.2 |
|
452 |
yygao |
1.5 |
// TagVertex Requirement
|
453 |
|
|
if ( isGoodTagVertex( *(splitVertexColl2->begin()), *(vertexColl->begin()), zsigncut_) ) {
|
454 |
yygao |
1.1 |
if(verbose_) {
|
455 |
yygao |
1.5 |
std::cout<<"splitTrackCollection1Handle->size() = " << int(splitTrackCollection1Handle->size()) << std::endl;
|
456 |
|
|
std::cout<<"splitTrackCollection2Handle->size() = " << int(splitTrackCollection2Handle->size()) << std::endl;
|
457 |
yygao |
1.1 |
}
|
458 |
yygao |
1.4 |
int nGoodTracks = 0;
|
459 |
yygao |
1.5 |
for(unsigned int i = 0; i < splitTrackCollection1Handle->size(); i++ ) {
|
460 |
|
|
edm::RefToBase<reco::Track> track(splitTrackCollection1Handle, i);
|
461 |
|
|
if(track->pt() > ptcut_) // && abs(track->eta()) < 2.5 )
|
462 |
yygao |
1.4 |
nGoodTracks++;
|
463 |
|
|
}
|
464 |
|
|
|
465 |
|
|
h_summary->Fill1d("denom_ntrack", nGoodTracks);
|
466 |
yygao |
1.5 |
// ProbeVertex Requirement
|
467 |
|
|
if ( isGoodProbeVertex( *splitVertexColl1->begin(), *vertexColl->begin(), zsigncut_) ) {
|
468 |
yygao |
1.4 |
h_summary->Fill1d("numer_ntrack", nGoodTracks);
|
469 |
yygao |
1.5 |
h_summary->Fill1d("avgWeight_orgvtx_eff",avgWeight(*vertexColl->begin()) );
|
470 |
|
|
h_summary->Fill1d("avgWeight_tagvtx_eff",avgWeight(*splitVertexColl2->begin()) );
|
471 |
|
|
h_summary->Fill1d("avgWeight_probevtx_eff",avgWeight(*splitVertexColl1->begin()) );
|
472 |
|
|
}
|
473 |
|
|
else
|
474 |
|
|
{
|
475 |
|
|
if (nGoodTracks >= 2) {
|
476 |
|
|
h_summary->Fill1d("avgWeight_orgvtx_ineff",avgWeight(*vertexColl->begin()) );
|
477 |
|
|
h_summary->Fill1d("avgWeight_tagvtx_ineff",avgWeight(*splitVertexColl2->begin()) );
|
478 |
|
|
h_summary->Fill1d("avgWeight_probevtx_ineff",avgWeight(*splitVertexColl1->begin()) );
|
479 |
|
|
|
480 |
|
|
std::cout<<"**********Inefficiency Found************** " <<std::endl;
|
481 |
|
|
std::cout<<"==== Original Vertex: " <<std::endl;
|
482 |
|
|
//printRecVtx(*vertexColl->begin());
|
483 |
|
|
printRecVtxs(vertexCollectionHandle);
|
484 |
|
|
std::cout<<"==== Probe Vertex: " <<std::endl;
|
485 |
|
|
//printRecVtx(*splitVertexColl1->begin());
|
486 |
|
|
printRecVtxs(splitVertexCollection1Handle);
|
487 |
|
|
std::cout<<"==== Tag Vertex: " <<std::endl;
|
488 |
|
|
//printRecVtx(*splitVertexColl2->begin());
|
489 |
|
|
printRecVtxs(splitVertexCollection2Handle);
|
490 |
|
|
std::cout<<"==== Probe cluster information "<< std::endl;
|
491 |
|
|
printCluster(clusters1, beamSpot);
|
492 |
|
|
std::cout<<"==== Tag cluster information "<< std::endl;
|
493 |
|
|
printCluster(clusters2, beamSpot);
|
494 |
|
|
std::cout<<"**********Inefficiency Found************** " <<std::endl;
|
495 |
|
|
}
|
496 |
|
|
}
|
497 |
yygao |
1.1 |
}
|
498 |
|
|
}
|
499 |
|
|
|
500 |
|
|
|
501 |
|
|
// ------------ method called once each job just before starting event loop ------------
|
502 |
|
|
void
|
503 |
|
|
PVEffAnalyzer::beginJob()
|
504 |
|
|
{
|
505 |
|
|
}
|
506 |
|
|
|
507 |
|
|
// ------------ method called once each job just after ending the event loop ------------
|
508 |
|
|
void
|
509 |
|
|
PVEffAnalyzer::endJob() {
|
510 |
|
|
edm::LogInfo("Analysis") << "[endJob] Analyzing PV info" << endl;
|
511 |
|
|
|
512 |
|
|
for (vector<int>::const_iterator it= datamode.begin(); it != datamode.end() ; ++it) {
|
513 |
|
|
string suffix;
|
514 |
|
|
edm::LogInfo("Debug")<<"datamode = "<< *it<<endl;
|
515 |
|
|
switch(*it) {
|
516 |
|
|
case 0: suffix = "";
|
517 |
|
|
break;
|
518 |
|
|
case 1: suffix = "_mct";
|
519 |
|
|
break;
|
520 |
|
|
}
|
521 |
yygao |
1.2 |
if (analyze_)
|
522 |
|
|
MakeEff(h_summary->ReadHisto1D(TString("numer_ntrack"+suffix)), h_summary->ReadHisto1D(TString("denom_ntrack"+suffix)), h_summary->ReadHisto1D(TString("eff_ntrack"+suffix)), false, 1);
|
523 |
yygao |
1.1 |
}
|
524 |
|
|
}
|
525 |
|
|
|
526 |
|
|
|
527 |
|
|
|
528 |
|
|
bool PVEffAnalyzer::isResonance(const HepMC::GenParticle * p){
|
529 |
|
|
double ctau=(pdt->particle( abs(p->pdg_id()) ))->lifetime();
|
530 |
|
|
edm::LogInfo("Debug") << "[isResonance] isResonance " << p->pdg_id() << " " << ctau << endl;
|
531 |
|
|
return ctau >0 && ctau <1e-6;
|
532 |
|
|
}
|
533 |
|
|
|
534 |
|
|
bool PVEffAnalyzer::isFinalstateParticle(const HepMC::GenParticle * p){
|
535 |
|
|
return ( !p->end_vertex() && p->status()==1 );
|
536 |
|
|
}
|
537 |
|
|
|
538 |
|
|
bool PVEffAnalyzer::isCharged(const HepMC::GenParticle * p){
|
539 |
|
|
const ParticleData * part = pdt->particle( p->pdg_id() );
|
540 |
|
|
if (part){
|
541 |
|
|
return part->charge()!=0;
|
542 |
|
|
}else{
|
543 |
|
|
// the new/improved particle table doesn't know anti-particles
|
544 |
|
|
return pdt->particle( -p->pdg_id() )!=0;
|
545 |
|
|
}
|
546 |
|
|
}
|
547 |
|
|
|
548 |
|
|
void PVEffAnalyzer::printSimVtxs(const Handle<SimVertexContainer> simVtxs){
|
549 |
|
|
int i=0;
|
550 |
|
|
for(SimVertexContainer::const_iterator vsim=simVtxs->begin();
|
551 |
|
|
vsim!=simVtxs->end(); ++vsim){
|
552 |
|
|
cout << i++ << ")" << scientific
|
553 |
|
|
<< " evtid=" << vsim->eventId().event()
|
554 |
|
|
<< " sim x=" << vsim->position().x()
|
555 |
|
|
<< " sim y=" << vsim->position().y()
|
556 |
|
|
<< " sim z=" << vsim->position().z()
|
557 |
|
|
<< " sim t=" << vsim->position().t()
|
558 |
|
|
<< " parent=" << vsim->parentIndex()
|
559 |
|
|
<< endl;
|
560 |
|
|
}
|
561 |
|
|
}
|
562 |
|
|
|
563 |
|
|
void PVEffAnalyzer::printRecVtxs(const Handle<reco::VertexCollection> recVtxs){
|
564 |
|
|
int ivtx=0;
|
565 |
|
|
for(reco::VertexCollection::const_iterator v=recVtxs->begin();
|
566 |
|
|
v!=recVtxs->end(); ++v){
|
567 |
|
|
cout << "Recvtx "<< std::setw(3) << std::setfill(' ')<<ivtx++
|
568 |
|
|
<< "#trk " << std::setw(3) << v->tracksSize()
|
569 |
|
|
<< " chi2 " << std::setw(4) << v->chi2()
|
570 |
|
|
<< " ndof " << std::setw(3) << v->ndof() << endl
|
571 |
|
|
<< " x " << std::setw(8) <<std::fixed << std::setprecision(4) << v->x()
|
572 |
|
|
<< " dx " << std::setw(8) << v->xError()<< endl
|
573 |
|
|
<< " y " << std::setw(8) << v->y()
|
574 |
|
|
<< " dy " << std::setw(8) << v->yError()<< endl
|
575 |
|
|
<< " z " << std::setw(8) << v->z()
|
576 |
|
|
<< " dz " << std::setw(8) << v->zError()
|
577 |
|
|
<< endl;
|
578 |
|
|
}
|
579 |
|
|
}
|
580 |
|
|
|
581 |
|
|
void PVEffAnalyzer::printRecVtx(const reco::Vertex & v){
|
582 |
|
|
|
583 |
|
|
cout << "#trk " << std::setw(3) << v.tracksSize()
|
584 |
|
|
<< " chi2 " << std::setw(4) << v.chi2()
|
585 |
|
|
<< " ndof " << std::setw(3) << v.ndof() << endl
|
586 |
|
|
<< " x " << std::setw(8) <<std::fixed << std::setprecision(4) << v.x()
|
587 |
|
|
<< " dx " << std::setw(8) << v.xError()<< endl
|
588 |
|
|
<< " y " << std::setw(8) << v.y()
|
589 |
|
|
<< " dy " << std::setw(8) << v.yError()<< endl
|
590 |
|
|
<< " z " << std::setw(8) << v.z()
|
591 |
|
|
<< " dz " << std::setw(8) << v.zError()
|
592 |
|
|
<< endl;
|
593 |
|
|
}
|
594 |
|
|
|
595 |
|
|
void PVEffAnalyzer::printSimTrks(const Handle<SimTrackContainer> simTrks){
|
596 |
|
|
cout << " simTrks type, (momentum), vertIndex, genpartIndex" << endl;
|
597 |
|
|
int i=1;
|
598 |
|
|
for(SimTrackContainer::const_iterator t=simTrks->begin();
|
599 |
|
|
t!=simTrks->end(); ++t){
|
600 |
|
|
//HepMC::GenParticle* gp=evtMC->GetEvent()->particle( (*t).genpartIndex() );
|
601 |
|
|
cout << i++ << ")"
|
602 |
|
|
<< (*t)
|
603 |
|
|
<< " index="
|
604 |
|
|
<< (*t).genpartIndex();
|
605 |
|
|
//if (gp) {
|
606 |
|
|
// HepMC::GenVertex *gv=gp->production_vertex();
|
607 |
|
|
// cout << " genvertex =" << (*gv);
|
608 |
|
|
//}
|
609 |
|
|
cout << endl;
|
610 |
|
|
}
|
611 |
|
|
}
|
612 |
|
|
|
613 |
|
|
// Select based on the highest SumPtSq Vertex
|
614 |
|
|
// Require at least 6 tracks with track weight > ~ 0.8
|
615 |
|
|
// vtx.ndof = 2*Sum(trackWeight) - 3
|
616 |
|
|
|
617 |
|
|
bool PVEffAnalyzer::isGoodSplitEvent( const reco::Vertex & org_vtx)
|
618 |
|
|
{
|
619 |
yygao |
1.5 |
if ( org_vtx.tracksSize() >= 6 && org_vtx.ndof() >= 6 && avgWeight(org_vtx) > 0.75 ) return true;
|
620 |
|
|
//if ( org_vtx.tracksSize() >= 6 && org_vtx.ndof() >= 6 ) return true;
|
621 |
|
|
else
|
622 |
|
|
return false;
|
623 |
yygao |
1.1 |
}
|
624 |
|
|
|
625 |
|
|
// Comparing the tag vertex with the unsplit vertex position in z
|
626 |
|
|
bool PVEffAnalyzer::isGoodTagVertex( const reco::Vertex & tag_vtx, const reco::Vertex & org_vtx, double & zsign_cut)
|
627 |
|
|
{
|
628 |
|
|
if(tag_vtx.isValid () && !tag_vtx.isFake()) {
|
629 |
|
|
float zsign = TMath::Abs(tag_vtx.z() - org_vtx.z() ) / TMath::Max(tag_vtx.zError(), org_vtx.zError() );
|
630 |
yygao |
1.5 |
//float zsign = TMath::Abs(tag_vtx.z() - org_vtx.z() ) / sqrt(pow(tag_vtx.zError(),2)+pow(org_vtx.zError(),2) );
|
631 |
yygao |
1.1 |
return (zsign < zsign_cut) ;
|
632 |
|
|
}
|
633 |
|
|
else
|
634 |
|
|
return false;
|
635 |
|
|
}
|
636 |
|
|
|
637 |
|
|
// Comparing the probe vertex with the unsplit vertex position in z
|
638 |
|
|
bool PVEffAnalyzer::isGoodProbeVertex( const reco::Vertex & probe_vtx, const reco::Vertex & org_vtx, double & zsign_cut)
|
639 |
|
|
{
|
640 |
|
|
if(probe_vtx.isValid () && !probe_vtx.isFake()) {
|
641 |
|
|
float zsign = TMath::Abs(probe_vtx.z() - org_vtx.z() ) / TMath::Max(probe_vtx.zError(), org_vtx.zError() );
|
642 |
yygao |
1.5 |
//float zsign = TMath::Abs(probe_vtx.z() - org_vtx.z() ) / sqrt(pow(probe_vtx.zError(),2)+pow(org_vtx.zError(),2) );
|
643 |
yygao |
1.1 |
return (zsign < zsign_cut) ;
|
644 |
|
|
}
|
645 |
|
|
else
|
646 |
|
|
return false;
|
647 |
|
|
}
|
648 |
|
|
|
649 |
|
|
|
650 |
|
|
bool PVEffAnalyzer::isAssoVertex(const PVEffAnalyzer::simPrimaryVertex & vsim,
|
651 |
|
|
const reco::Vertex & vrec, double & zsign_cut)
|
652 |
|
|
{
|
653 |
|
|
if(vrec.isValid() && !vrec.isFake() )
|
654 |
|
|
return (TMath::Abs(vsim.z-vrec.z()) < zsign_cut * vrec.zError() );
|
655 |
|
|
else
|
656 |
|
|
return false;
|
657 |
|
|
}
|
658 |
|
|
|
659 |
|
|
|
660 |
yygao |
1.5 |
|
661 |
yygao |
1.1 |
void PVEffAnalyzer::MakeEff(TH1D* numer, TH1D* denom, TH1D* eff, //TGraphAsymmErrors* & gr_eff,
|
662 |
|
|
const bool rebin, const Float_t n ) {
|
663 |
|
|
eff->Divide( numer, denom, 1,1,"B" );
|
664 |
|
|
|
665 |
|
|
if( rebin ) {
|
666 |
|
|
std::vector<Double_t> binEdges;
|
667 |
|
|
Int_t bin = denom->GetNbinsX() + 1;
|
668 |
|
|
binEdges.push_back(denom->GetBinLowEdge(bin));
|
669 |
|
|
Float_t nEntries = 0;
|
670 |
|
|
while (bin > 1) {
|
671 |
|
|
bin --;
|
672 |
|
|
nEntries = denom->GetBinContent(bin);
|
673 |
|
|
while (nEntries < n && bin > 1) {
|
674 |
|
|
bin --;
|
675 |
|
|
nEntries += denom->GetBinContent(bin);
|
676 |
|
|
}
|
677 |
|
|
binEdges.push_back(denom->GetBinLowEdge(bin));
|
678 |
|
|
}
|
679 |
|
|
|
680 |
|
|
Double_t *array = new Double_t[binEdges.size()];
|
681 |
|
|
|
682 |
|
|
Int_t j = 0;
|
683 |
|
|
for (Int_t i = binEdges.size(); i > 0; --i) {
|
684 |
|
|
array[j] = binEdges[i - 1];
|
685 |
|
|
++j;
|
686 |
|
|
}
|
687 |
|
|
}
|
688 |
|
|
|
689 |
|
|
for (int i = 1; i<eff->GetNbinsX()+1; i++) {
|
690 |
|
|
if(numer->GetBinContent(i) == 0 || denom->GetBinContent(i) == 0 ) continue;
|
691 |
yygao |
1.3 |
float error = sqrt(eff->GetBinContent(i)*(1-eff->GetBinContent(i))/denom->GetBinContent(i)); // Binominal-Error
|
692 |
yygao |
1.1 |
eff->SetBinError(i, error);
|
693 |
|
|
}
|
694 |
|
|
}
|
695 |
yygao |
1.4 |
|
696 |
yygao |
1.5 |
bool PVEffAnalyzer::goodCluster(std::vector< std::vector<reco::TransientTrack> > & clusters)
|
697 |
|
|
{
|
698 |
|
|
if( clusters.size() == 1 && ( clusters.begin()->size() > 1 )) return true;
|
699 |
|
|
else
|
700 |
|
|
return false;
|
701 |
|
|
}
|
702 |
|
|
|
703 |
|
|
void PVEffAnalyzer::printCluster(std::vector< std::vector<reco::TransientTrack> > & clusters, const Point & beamSpot)
|
704 |
|
|
{
|
705 |
|
|
//std::cout<<"theTrackClusterizer.zSeparation()"<<theTrackClusterizer.zSeparation()<<std::endl;
|
706 |
|
|
cout << "#clusters " << std::setw(3) << clusters.size() << endl;
|
707 |
|
|
int idx_cluster = 0;
|
708 |
|
|
for (vector< vector<reco::TransientTrack> >::const_iterator iclus = clusters.begin(); iclus != clusters.end(); iclus++, idx_cluster++) {
|
709 |
|
|
if((*iclus).size() == 0) continue;
|
710 |
|
|
cout << "=== Cluster # " << std::setw(3) << idx_cluster << std::setw(3) << " nTracks" << std::setw(3) << (*iclus).size() <<endl;
|
711 |
|
|
for (int idx_track = 0; idx_track < (*iclus).size() ; idx_track++) {
|
712 |
|
|
reco::Track t = (*iclus).at(idx_track).track();
|
713 |
|
|
cout << "== Trk # " << std::setw(4) << idx_track
|
714 |
|
|
<< " nHit " << std::setw(4) << t.hitPattern().numberOfValidHits() <<endl
|
715 |
|
|
<< " pt " << std::setw(9) << std::fixed << std::setprecision(4) << t.pt()
|
716 |
|
|
<< " dpt " << std::setw(9) << t.ptError() << endl
|
717 |
|
|
<< " dxy(BS) " << std::setw(9) << t.dxy(beamSpot)
|
718 |
|
|
<< " sigmaDxy " << std::setw(9) << t.dxyError() << endl
|
719 |
|
|
<< " dz(BS) " << std::setw(9) << t.dz(beamSpot)
|
720 |
|
|
<< " sigmaDz " << std::setw(9) << t.dzError() << endl
|
721 |
|
|
<< " xPCA " << std::setw(9) << t.vertex().x()
|
722 |
|
|
<< " yPCA " << std::setw(9) << t.vertex().y()
|
723 |
|
|
<< " zPCA " << std::setw(9) << t.vertex().z()
|
724 |
|
|
<< endl;
|
725 |
|
|
}
|
726 |
|
|
}
|
727 |
|
|
|
728 |
|
|
}
|
729 |
|
|
double PVEffAnalyzer::avgWeight(const reco::Vertex & vtx)
|
730 |
|
|
{
|
731 |
|
|
if(!vtx.isValid() || vtx.isFake() )
|
732 |
|
|
return 0;
|
733 |
|
|
else
|
734 |
|
|
return (vtx.ndof()+3.0)/2.0/double(vtx.tracksSize());
|
735 |
|
|
}
|