Main Page
Namespaces
Classes
Package Documentation
RecoVertex
KalmanVertexFit
plugins
KVFTest.cc
Go to the documentation of this file.
1
#include "
RecoVertex/KalmanVertexFit/plugins/KVFTest.h
"
2
3
#include "
DataFormats/VertexReco/interface/Vertex.h
"
4
#include "
DataFormats/VertexReco/interface/VertexFwd.h
"
5
#include "
DataFormats/TrackReco/interface/Track.h
"
6
#include "
DataFormats/Common/interface/Handle.h
"
7
#include "
FWCore/Framework/interface/MakerMacros.h
"
8
#include "
FWCore/MessageLogger/interface/MessageLogger.h
"
9
#include "
FWCore/Framework/interface/ESHandle.h
"
10
11
#include "
TrackingTools/TransientTrack/interface/TransientTrackBuilder.h
"
12
#include "
TrackingTools/Records/interface/TransientTrackRecord.h
"
13
#include "
TrackingTools/TransientTrack/interface/TransientTrack.h
"
14
#include "
RecoVertex/VertexPrimitives/interface/TransientVertex.h
"
15
#include "
RecoVertex/KalmanVertexFit/interface/KalmanVertexFitter.h
"
16
#include "
MagneticField/Engine/interface/MagneticField.h
"
17
#include "
MagneticField/Records/interface/IdealMagneticFieldRecord.h
"
18
19
#include <iostream>
20
21
using namespace
reco
;
22
using namespace
edm
;
23
using namespace
std
;
24
25
KVFTest::KVFTest
(
const
edm::ParameterSet
& iConfig)
26
: theConfig(iConfig)
27
{
28
token_tracks
= consumes<TrackCollection>(iConfig.
getParameter
<
string
>(
"TrackLabel"
));
29
outputFile_
= iConfig.
getUntrackedParameter
<
std::string
>(
"outputFile"
);
30
kvfPSet
= iConfig.
getParameter
<
edm::ParameterSet
>(
"KVFParameters"
);
31
rootFile_
= TFile::Open(
outputFile_
.c_str(),
"RECREATE"
);
32
edm::LogInfo
(
"RecoVertex/KVFTest"
)
33
<<
"Initializing KVF TEST analyser - Output file: "
<<
outputFile_
<<
"\n"
;
34
35
token_TrackTruth
= consumes<TrackingParticleCollection>(
edm::InputTag
(
"trackingtruth"
,
"TrackTruth"
));
36
token_VertexTruth
= consumes<TrackingVertexCollection>(
edm::InputTag
(
"trackingtruth"
,
"VertexTruth"
));
37
token_associatorForParamAtPca
= consumes<reco::TrackToTrackingParticleAssociator>(
edm::InputTag
(
"trackAssociatorByChi2"
));
38
39
}
40
41
42
KVFTest::~KVFTest
() {
43
delete
rootFile_
;
44
}
45
46
void
KVFTest::beginJob
(){
47
}
48
49
50
void
KVFTest::endJob
() {
51
}
52
53
//
54
// member functions
55
//
56
57
void
58
KVFTest::analyze
(
const
edm::Event
&
iEvent
,
const
edm::EventSetup
& iSetup)
59
{
60
edm::Handle<reco::TrackToTrackingParticleAssociator>
associatorForParamAtPca;
61
iEvent.
getByToken
(
token_associatorForParamAtPca
,associatorForParamAtPca);
62
63
if
(not
tree
) {
64
edm::ESHandle<MagneticField>
magField;
65
iSetup.
get
<
IdealMagneticFieldRecord
>().
get
(magField);
66
tree.reset(
new
SimpleVertexTree
(
"VertexFitter"
, magField.
product
()) );
67
}
68
69
70
71
edm::LogInfo
(
"RecoVertex/KVFTest"
)
72
<<
"Reconstructing event number: "
<< iEvent.
id
() <<
"\n"
;
73
74
// get RECO tracks from the event
75
// `tks` can be used as a ptr to a reco::TrackCollection
76
edm::Handle<edm::View<reco::Track>
> tks;
77
iEvent.
getByToken
(
token_tracks
, tks);
78
if
(!tks.
isValid
()) {
79
edm::LogInfo
(
"RecoVertex/KVFTest"
)
80
<<
"Exception during event number: "
<< iEvent.
id
()
81
<<
"\n"
;
82
}
else
{
83
edm::LogInfo
(
"RecoVertex/KVFTest"
)
84
<<
"Found: "
<< (*tks).size() <<
" reconstructed tracks"
<<
"\n"
;
85
std::cout
<<
"got "
<< (*tks).size() <<
" tracks "
<< std::endl;
86
87
// Transform Track to TransientTrack
88
89
//get the builder:
90
edm::ESHandle<TransientTrackBuilder>
theB;
91
iSetup.
get
<
TransientTrackRecord
>().
get
(
"TransientTrackBuilder"
,theB);
92
//do the conversion:
93
std::vector<TransientTrack> t_tks = (*theB).build(tks);
94
95
edm::LogInfo
(
"RecoVertex/KVFTest"
)
96
<<
"Found: "
<< t_tks.size() <<
" reconstructed tracks"
<<
"\n"
;
97
98
// Call the KalmanVertexFitter if more than 1 track
99
if
(t_tks.size() > 1) {
100
// KalmanVertexFitter kvf(kvfPSet);
101
KalmanVertexFitter
kvf(
true
);
102
TransientVertex
tv = kvf.
vertex
(t_tks);
103
104
std::cout
<<
"Position: "
<<
Vertex::Point
(tv.
position
()) <<
"\n"
;
105
106
// For the analysis: compare to your SimVertex
107
TrackingVertex
sv =
getSimVertex
(iEvent);
108
edm::Handle<TrackingParticleCollection>
TPCollectionH ;
109
iEvent.
getByToken
(
token_TrackTruth
, TPCollectionH);
110
if
(!TPCollectionH.
isValid
()) {
111
edm::LogInfo
(
"RecoVertex/KVFTest"
)
112
<<
"Exception during event number: "
<< iEvent.
id
()
113
<<
"\n"
;
114
}
else
{
115
const
TrackingParticleCollection
tPC = *(TPCollectionH.
product
());
116
reco::RecoToSimCollection
recSimColl=associatorForParamAtPca->
associateRecoToSim
(tks,
117
TPCollectionH);
118
tree->fill(tv, &sv, &recSimColl);
119
}
120
}
121
}
122
}
123
124
//Returns the first vertex in the list.
125
126
TrackingVertex
KVFTest::getSimVertex
(
const
edm::Event
&
iEvent
)
const
127
{
128
// get the simulated vertices
129
edm::Handle<TrackingVertexCollection>
TVCollectionH ;
130
iEvent.
getByToken
(
token_VertexTruth
,TVCollectionH);
131
const
TrackingVertexCollection
tPC = *(TVCollectionH.
product
());
132
133
// Handle<edm::SimVertexContainer> simVtcs;
134
// iEvent.getByLabel("g4SimHits", simVtcs);
135
// std::cout << "SimVertex " << simVtcs->size() << std::endl;
136
// for(edm::SimVertexContainer::const_iterator v=simVtcs->begin();
137
// v!=simVtcs->end(); ++v){
138
// std::cout << "simvtx "
139
// << v->position().x() << " "
140
// << v->position().y() << " "
141
// << v->position().z() << " "
142
// << v->parentIndex() << " "
143
// << v->noParent() << " "
144
// << std::endl;
145
// }
146
return
*(tPC.begin());
147
}
148
DEFINE_FWK_MODULE
(
KVFTest
);
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
KVFTest::outputFile_
std::string outputFile_
Definition:
KVFTest.h:56
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
TransientTrackBuilder.h
MessageLogger.h
KVFTest
Definition:
KVFTest.h:37
TrackingVertex
Definition:
TrackingVertex.h:22
TrackingParticleCollection
std::vector< TrackingParticle > TrackingParticleCollection
Definition:
HiEvtPlaneFlatProducer.cc:111
KVFTest::getSimVertex
TrackingVertex getSimVertex(const edm::Event &iEvent) const
Definition:
KVFTest.cc:126
edm::Event::getByToken
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition:
Event.h:519
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition:
AlCaHLTBitMon_QueryRunRegistry.py:255
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition:
MakerMacros.h:17
TransientTrack.h
KalmanVertexFitter
Definition:
KalmanVertexFitter.h:22
KVFTest::token_tracks
edm::EDGetTokenT< reco::TrackCollection > token_tracks
Definition:
KVFTest.h:57
MakerMacros.h
edm::Handle< reco::TrackToTrackingParticleAssociator >
IdealMagneticFieldRecord
Definition:
IdealMagneticFieldRecord.h:11
std
Definition:
JetResolutionObject.h:76
VertexFwd.h
pftools::Point
std::pair< double, double > Point
Definition:
CaloEllipse.h:18
KVFTest::token_associatorForParamAtPca
edm::EDGetTokenT< reco::TrackToTrackingParticleAssociator > token_associatorForParamAtPca
Definition:
KVFTest.h:60
MagneticField.h
SimpleVertexTree
Definition:
SimpleVertexTree.h:41
KVFTest::token_TrackTruth
edm::EDGetTokenT< TrackingParticleCollection > token_TrackTruth
Definition:
KVFTest.h:58
iEvent
int iEvent
Definition:
GenABIO.cc:230
TransientVertex.h
edm::ESHandle< MagneticField >
TransientVertex::position
GlobalPoint position() const
Definition:
TransientVertex.h:136
ESHandle.h
KVFTest::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition:
KVFTest.cc:58
edm::AssociationMap< edm::OneToManyWithQualityGeneric< edm::View< reco::Track >, TrackingParticleCollection, double > >
KalmanVertexFitter::vertex
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const override
Definition:
KalmanVertexFitter.h:59
KVFTest::rootFile_
TFile * rootFile_
Definition:
KVFTest.h:54
edm::EventSetup
Definition:
EventSetup.h:48
edm::HandleBase::isValid
bool isValid() const
Definition:
HandleBase.h:74
KVFTest::~KVFTest
~KVFTest() override
Definition:
KVFTest.cc:42
Vertex.h
edm::LogInfo
Definition:
MessageLogger.h:216
edm::Handle::product
T const * product() const
Definition:
Handle.h:81
KVFTest::token_VertexTruth
edm::EDGetTokenT< TrackingVertexCollection > token_VertexTruth
Definition:
KVFTest.h:59
reco::TrackToTrackingParticleAssociator::associateRecoToSim
reco::RecoToSimCollection associateRecoToSim(const edm::Handle< edm::View< reco::Track > > &tCH, const edm::Handle< TrackingParticleCollection > &tPCH) const
compare reco to sim the handle of reco::Track and TrackingParticle collections
Definition:
TrackToTrackingParticleAssociator.h:65
KalmanVertexFitter.h
KVFTest::beginJob
void beginJob() override
Definition:
KVFTest.cc:46
TrackingVertexCollection
std::vector< TrackingVertex > TrackingVertexCollection
Definition:
TrackingVertexContainer.h:8
edm::EventSetup::get
const T & get() const
Definition:
EventSetup.h:59
TransientVertex
Definition:
TransientVertex.h:18
KVFTest::KVFTest
KVFTest(const edm::ParameterSet &)
Definition:
KVFTest.cc:25
TransientTrackRecord
Definition:
TransientTrackRecord.h:12
edm::EventBase::id
edm::EventID id() const
Definition:
EventBase.h:60
TransientTrackRecord.h
reco
fixed size matrix
Definition:
AlignmentAlgorithmBase.h:43
edm
HLT enums.
Definition:
AlignableModifier.h:17
edm::InputTag
Definition:
InputTag.h:15
KVFTest::endJob
void endJob() override
Definition:
KVFTest.cc:50
KVFTest::kvfPSet
edm::ParameterSet kvfPSet
Definition:
KVFTest.h:52
edm::ParameterSet
Definition:
ParameterSet.h:36
gather_cfg.cout
cout
Definition:
gather_cfg.py:145
tree
Definition:
tree.py:1
edm::Event
Definition:
Event.h:69
Track.h
edm::ESHandle::product
T const * product() const
Definition:
ESHandle.h:86
IdealMagneticFieldRecord.h
Handle.h
KVFTest.h
Generated for CMSSW Reference Manual by
1.8.11