1 |
dgele |
1.1 |
//
|
2 |
|
|
// $Id: LeptonVertexSignificance.cc,v 1.2 2008/02/28 14:54:25 llista Exp $
|
3 |
|
|
//
|
4 |
|
|
|
5 |
|
|
#include "PhysicsTools/PatUtils/interface/LeptonVertexSignificance.h"
|
6 |
|
|
#include "FWCore/MessageLogger/interface/MessageLogger.h"
|
7 |
|
|
#include "FWCore/Utilities/interface/Exception.h"
|
8 |
|
|
#include "FWCore/Framework/interface/ESHandle.h"
|
9 |
|
|
#include "DataFormats/Common/interface/Handle.h"
|
10 |
|
|
#include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
|
11 |
|
|
#include "TrackingTools/Records/interface/TransientTrackRecord.h"
|
12 |
|
|
#include "TrackingTools/TransientTrack/interface/TransientTrack.h"
|
13 |
|
|
#include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h"
|
14 |
|
|
#include "DataFormats/VertexReco/interface/Vertex.h"
|
15 |
|
|
#include "DataFormats/VertexReco/interface/VertexFwd.h"
|
16 |
|
|
#include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
|
17 |
|
|
#include "FWCore/Framework/interface/EventSetup.h"
|
18 |
|
|
#include "FWCore/Framework/interface/Event.h"
|
19 |
|
|
#include "DataFormats/PatCandidates/interface/Electron.h"
|
20 |
|
|
#include "DataFormats/PatCandidates/interface/Muon.h"
|
21 |
|
|
|
22 |
|
|
using namespace pat;
|
23 |
|
|
|
24 |
|
|
// constructor
|
25 |
|
|
LeptonVertexSignificance::LeptonVertexSignificance(const edm::EventSetup & iSetup) {
|
26 |
|
|
// instantiate the transient-track builder
|
27 |
|
|
edm::ESHandle<TransientTrackBuilder> builder;
|
28 |
|
|
iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder", builder);
|
29 |
|
|
theTrackBuilder_ = new TransientTrackBuilder(*builder.product());
|
30 |
|
|
}
|
31 |
|
|
|
32 |
|
|
// destructor
|
33 |
|
|
LeptonVertexSignificance::~LeptonVertexSignificance() {
|
34 |
|
|
delete theTrackBuilder_;
|
35 |
|
|
}
|
36 |
|
|
|
37 |
|
|
// calculate the TrackIsoPt for the lepton object
|
38 |
|
|
float LeptonVertexSignificance::calculate(const Electron & theElectron, const edm::Event & iEvent) {
|
39 |
|
|
return this->calculate(*theElectron.gsfTrack(), iEvent);
|
40 |
|
|
}
|
41 |
|
|
|
42 |
|
|
float LeptonVertexSignificance::calculate(const Muon & theMuon, const edm::Event & iEvent) {
|
43 |
|
|
return this->calculate(*theMuon.track(), iEvent);
|
44 |
|
|
}
|
45 |
|
|
|
46 |
|
|
// calculate the TrackIsoPt for the lepton's track
|
47 |
|
|
float LeptonVertexSignificance::calculate(const reco::Track & theTrack, const edm::Event & iEvent) {
|
48 |
|
|
// FIXME: think more about how to handle events without vertices
|
49 |
|
|
// lepton LR calculation should have nothing to do with event selection
|
50 |
|
|
edm::Handle<reco::VertexCollection> vertexHandle;
|
51 |
|
|
iEvent.getByLabel("offlinePrimaryVerticesFromCTFTracks", vertexHandle);
|
52 |
|
|
if (vertexHandle.product()->size() == 0) return 0;
|
53 |
|
|
reco::Vertex theVertex = vertexHandle.product()->front();
|
54 |
|
|
// calculate the track-vertex association significance
|
55 |
|
|
reco::TransientTrack theTrTrack = theTrackBuilder_->build(&theTrack);
|
56 |
|
|
GlobalPoint theVertexPoint(theVertex.position().x(), theVertex.position().y(), theVertex.position().z());
|
57 |
|
|
FreeTrajectoryState theLeptonNearVertex = theTrTrack.trajectoryStateClosestToPoint(theVertexPoint).theState();
|
58 |
|
|
return fabs(theVertex.position().z() - theLeptonNearVertex.position().z())
|
59 |
|
|
/ sqrt(pow(theVertex.zError(), 2) + theLeptonNearVertex.cartesianError().position().czz());
|
60 |
|
|
}
|
61 |
|
|
|