CMS 3D CMS Logo

KVFTest.cc
Go to the documentation of this file.
2 
10 
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) : theConfig(iConfig) {
26  token_tracks = consumes<TrackCollection>(iConfig.getParameter<string>("TrackLabel"));
27  outputFile_ = iConfig.getUntrackedParameter<std::string>("outputFile");
28  kvfPSet = iConfig.getParameter<edm::ParameterSet>("KVFParameters");
29  rootFile_ = TFile::Open(outputFile_.c_str(), "RECREATE");
30  edm::LogInfo("RecoVertex/KVFTest") << "Initializing KVF TEST analyser - Output file: " << outputFile_ << "\n";
31 
32  token_TrackTruth = consumes<TrackingParticleCollection>(edm::InputTag("trackingtruth", "TrackTruth"));
33  token_VertexTruth = consumes<TrackingVertexCollection>(edm::InputTag("trackingtruth", "VertexTruth"));
35  consumes<reco::TrackToTrackingParticleAssociator>(edm::InputTag("trackAssociatorByChi2"));
36 }
37 
39 
41 
42 void KVFTest::endJob() {}
43 
44 //
45 // member functions
46 //
47 
48 void KVFTest::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
50  iEvent.getByToken(token_associatorForParamAtPca, associatorForParamAtPca);
51 
52  if (not tree) {
54  iSetup.get<IdealMagneticFieldRecord>().get(magField);
55  tree.reset(new SimpleVertexTree("VertexFitter", magField.product()));
56  }
57 
58  edm::LogInfo("RecoVertex/KVFTest") << "Reconstructing event number: " << iEvent.id() << "\n";
59 
60  // get RECO tracks from the event
61  // `tks` can be used as a ptr to a reco::TrackCollection
63  iEvent.getByToken(token_tracks, tks);
64  if (!tks.isValid()) {
65  edm::LogInfo("RecoVertex/KVFTest") << "Exception during event number: " << iEvent.id() << "\n";
66  } else {
67  edm::LogInfo("RecoVertex/KVFTest") << "Found: " << (*tks).size() << " reconstructed tracks"
68  << "\n";
69  std::cout << "got " << (*tks).size() << " tracks " << std::endl;
70 
71  // Transform Track to TransientTrack
72 
73  //get the builder:
75  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder", theB);
76  //do the conversion:
77  std::vector<TransientTrack> t_tks = (*theB).build(tks);
78 
79  edm::LogInfo("RecoVertex/KVFTest") << "Found: " << t_tks.size() << " reconstructed tracks"
80  << "\n";
81 
82  // Call the KalmanVertexFitter if more than 1 track
83  if (t_tks.size() > 1) {
84  // KalmanVertexFitter kvf(kvfPSet);
85  KalmanVertexFitter kvf(true);
86  TransientVertex tv = kvf.vertex(t_tks);
87 
88  std::cout << "Position: " << Vertex::Point(tv.position()) << "\n";
89 
90  // For the analysis: compare to your SimVertex
93  iEvent.getByToken(token_TrackTruth, TPCollectionH);
94  if (!TPCollectionH.isValid()) {
95  edm::LogInfo("RecoVertex/KVFTest") << "Exception during event number: " << iEvent.id() << "\n";
96  } else {
97  const TrackingParticleCollection tPC = *(TPCollectionH.product());
98  reco::RecoToSimCollection recSimColl = associatorForParamAtPca->associateRecoToSim(tks, TPCollectionH);
99  tree->fill(tv, &sv, &recSimColl);
100  }
101  }
102  }
103 }
104 
105 //Returns the first vertex in the list.
106 
108  // get the simulated vertices
110  iEvent.getByToken(token_VertexTruth, TVCollectionH);
111  const TrackingVertexCollection tPC = *(TVCollectionH.product());
112 
113  // Handle<edm::SimVertexContainer> simVtcs;
114  // iEvent.getByLabel("g4SimHits", simVtcs);
115  // std::cout << "SimVertex " << simVtcs->size() << std::endl;
116  // for(edm::SimVertexContainer::const_iterator v=simVtcs->begin();
117  // v!=simVtcs->end(); ++v){
118  // std::cout << "simvtx "
119  // << v->position().x() << " "
120  // << v->position().y() << " "
121  // << v->position().z() << " "
122  // << v->parentIndex() << " "
123  // << v->noParent() << " "
124  // << std::endl;
125  // }
126  return *(tPC.begin());
127 }
edm::ESHandle::product
T const * product() const
Definition: ESHandle.h:86
KalmanVertexFitter::vertex
CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const override
Definition: KalmanVertexFitter.h:49
Handle.h
KVFTest::KVFTest
KVFTest(const edm::ParameterSet &)
Definition: KVFTest.cc:25
MessageLogger.h
KVFTest::rootFile_
TFile * rootFile_
Definition: KVFTest.h:52
edm::Handle::product
T const * product() const
Definition: Handle.h:70
KalmanVertexFitter.h
ESHandle.h
edm
HLT enums.
Definition: AlignableModifier.h:19
tree
Definition: tree.py:1
reco::TrackToTrackingParticleAssociator::associateRecoToSim
reco::RecoToSimCollection associateRecoToSim(const edm::Handle< edm::View< reco::Track >> &tCH, const edm::Handle< TrackingParticleCollection > &tPCH) const
Definition: TrackToTrackingParticleAssociator.h:64
gather_cfg.cout
cout
Definition: gather_cfg.py:144
edm::LogInfo
Definition: MessageLogger.h:254
TransientVertex::position
GlobalPoint position() const
Definition: TransientVertex.h:169
TrackingVertexCollection
std::vector< TrackingVertex > TrackingVertexCollection
Definition: TrackingVertexContainer.h:8
KVFTest::~KVFTest
~KVFTest() override
Definition: KVFTest.cc:38
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
TransientTrack.h
edm::Handle< reco::TrackToTrackingParticleAssociator >
IdealMagneticFieldRecord
Definition: IdealMagneticFieldRecord.h:11
KVFTest::token_associatorForParamAtPca
edm::EDGetTokenT< reco::TrackToTrackingParticleAssociator > token_associatorForParamAtPca
Definition: KVFTest.h:58
MakerMacros.h
Point
math::XYZPoint Point
Definition: TrackerDpgAnalysis.cc:107
Track.h
edm::EventSetup::get
T get() const
Definition: EventSetup.h:73
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
KVFTest::token_tracks
edm::EDGetTokenT< reco::TrackCollection > token_tracks
Definition: KVFTest.h:55
pfDeepBoostedJetPreprocessParams_cfi.sv
sv
Definition: pfDeepBoostedJetPreprocessParams_cfi.py:226
TransientTrackRecord
Definition: TransientTrackRecord.h:11
IdealMagneticFieldRecord.h
edm::ESHandle< MagneticField >
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
Vertex.h
KVFTest.h
KVFTest::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: KVFTest.cc:48
KVFTest::beginJob
void beginJob() override
Definition: KVFTest.cc:40
TransientTrackBuilder.h
edm::ParameterSet
Definition: ParameterSet.h:36
KVFTest::endJob
void endJob() override
Definition: KVFTest.cc:42
edm::AssociationMap< edm::OneToManyWithQualityGeneric< edm::View< reco::Track >, TrackingParticleCollection, double > >
TrackingVertex
Definition: TrackingVertex.h:22
KVFTest
Definition: KVFTest.h:36
iEvent
int iEvent
Definition: GenABIO.cc:224
KVFTest::kvfPSet
edm::ParameterSet kvfPSet
Definition: KVFTest.h:50
MagneticField.h
TransientVertex
Definition: TransientVertex.h:18
edm::EventSetup
Definition: EventSetup.h:57
TransientTrackRecord.h
get
#define get
KVFTest::token_TrackTruth
edm::EDGetTokenT< TrackingParticleCollection > token_TrackTruth
Definition: KVFTest.h:56
VertexFwd.h
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
TransientVertex.h
std
Definition: JetResolutionObject.h:76
KVFTest::getSimVertex
TrackingVertex getSimVertex(const edm::Event &iEvent) const
Definition: KVFTest.cc:107
TrackingParticleCollection
std::vector< TrackingParticle > TrackingParticleCollection
Definition: TrackingParticleFwd.h:8
SimpleVertexTree
Definition: SimpleVertexTree.h:41
KVFTest::token_VertexTruth
edm::EDGetTokenT< TrackingVertexCollection > token_VertexTruth
Definition: KVFTest.h:57
edm::HandleBase::isValid
bool isValid() const
Definition: HandleBase.h:70
edm::Event
Definition: Event.h:73
KVFTest::outputFile_
std::string outputFile_
Definition: KVFTest.h:54
edm::InputTag
Definition: InputTag.h:15
KalmanVertexFitter
Definition: KalmanVertexFitter.h:22