CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
KVFTest.cc
Go to the documentation of this file.
2 
11 
18 
19 #include <iostream>
20 
21 using namespace reco;
22 using namespace edm;
23 using namespace std;
24 
26  : theConfig(iConfig), associatorForParamAtPca(0), tree(0)
27 {
28  trackLabel_ = iConfig.getParameter<std::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 
36 
38  delete rootFile_;
39 }
40 
42 }
43 
44 
46  delete tree;
47 }
48 
49 //
50 // member functions
51 //
52 
53 void
55 {
56  if ( associatorForParamAtPca==0 ) {
57  edm::ESHandle<TrackAssociatorBase> theAssociatorForParamAtPca;
58  iSetup.get<TrackAssociatorRecord>().get("TrackAssociatorByChi2",theAssociatorForParamAtPca);
59  associatorForParamAtPca = (TrackAssociatorByChi2 *) theAssociatorForParamAtPca.product();
60 
61  tree = new SimpleVertexTree("VertexFitter", associatorForParamAtPca);
62  }
63 
64 
65 
66  edm::LogInfo("RecoVertex/KVFTest")
67  << "Reconstructing event number: " << iEvent.id() << "\n";
68 
69  // get RECO tracks from the event
70  // `tks` can be used as a ptr to a reco::TrackCollection
72  iEvent.getByLabel(trackLabel_, tks);
73  if (!tks.isValid()) {
74  edm::LogInfo("RecoVertex/KVFTest")
75  << "Exception during event number: " << iEvent.id()
76  << "\n";
77  } else {
78  edm::LogInfo("RecoVertex/KVFTest")
79  << "Found: " << (*tks).size() << " reconstructed tracks" << "\n";
80  std::cout << "got " << (*tks).size() << " tracks " << std::endl;
81 
82  // Transform Track to TransientTrack
83 
84  //get the builder:
86  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
87  //do the conversion:
88  std::vector<TransientTrack> t_tks = (*theB).build(tks);
89 
90  edm::LogInfo("RecoVertex/KVFTest")
91  << "Found: " << t_tks.size() << " reconstructed tracks" << "\n";
92 
93  // Call the KalmanVertexFitter if more than 1 track
94  if (t_tks.size() > 1) {
95  // KalmanVertexFitter kvf(kvfPSet);
96  KalmanVertexFitter kvf(true);
97  TransientVertex tv = kvf.vertex(t_tks);
98 
99  std::cout << "Position: " << Vertex::Point(tv.position()) << "\n";
100 
101  // For the analysis: compare to your SimVertex
102  TrackingVertex sv = getSimVertex(iEvent);
104  iEvent.getByLabel("trackingtruth","TrackTruth",TPCollectionH);
105  if (!TPCollectionH.isValid()) {
106  edm::LogInfo("RecoVertex/KVFTest")
107  << "Exception during event number: " << iEvent.id()
108  << "\n";
109  } else {
110  const TrackingParticleCollection tPC = *(TPCollectionH.product());
112  TPCollectionH,
113  &iEvent);
114  tree->fill(tv, &sv, &recSimColl);
115  }
116  }
117  }
118 }
119 
120 //Returns the first vertex in the list.
121 
123 {
124  // get the simulated vertices
126  iEvent.getByLabel("trackingtruth","VertexTruth",TVCollectionH);
127  const TrackingVertexCollection tPC = *(TVCollectionH.product());
128 
129 // Handle<edm::SimVertexContainer> simVtcs;
130 // iEvent.getByLabel("g4SimHits", simVtcs);
131 // std::cout << "SimVertex " << simVtcs->size() << std::endl;
132 // for(edm::SimVertexContainer::const_iterator v=simVtcs->begin();
133 // v!=simVtcs->end(); ++v){
134 // std::cout << "simvtx "
135 // << v->position().x() << " "
136 // << v->position().y() << " "
137 // << v->position().z() << " "
138 // << v->parentIndex() << " "
139 // << v->noParent() << " "
140 // << std::endl;
141 // }
142  return *(tPC.begin());
143 }
TrackAssociatorByChi2 * associatorForParamAtPca
Definition: KVFTest.h:58
T getParameter(std::string const &) const
std::string outputFile_
Definition: KVFTest.h:62
T getUntrackedParameter(std::string const &, T const &) const
std::vector< TrackingParticle > TrackingParticleCollection
TrackingVertex getSimVertex(const edm::Event &iEvent) const
Definition: KVFTest.cc:122
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
SimpleVertexTree * tree
Definition: KVFTest.h:59
virtual CachingVertex< 5 > vertex(const std::vector< reco::TransientTrack > &tracks) const
void fill(const TransientVertex &recv, const TrackingVertex *simv=0, reco::RecoToSimCollection *recSimColl=0, const float &time=0.)
std::pair< double, double > Point
Definition: CaloEllipse.h:18
virtual void analyze(const edm::Event &, const edm::EventSetup &)
Definition: KVFTest.cc:54
virtual void endJob()
Definition: KVFTest.cc:45
int iEvent
Definition: GenABIO.cc:243
GlobalPoint position() const
std::string trackLabel_
Definition: KVFTest.h:63
TFile * rootFile_
Definition: KVFTest.h:60
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
double associateRecoToSim(reco::TrackCollection::const_iterator, TrackingParticleCollection::const_iterator, const reco::BeamSpot &) const
compare reco::TrackCollection and TrackingParticleCollection iterators: returns the chi2 ...
virtual void beginJob()
Definition: KVFTest.cc:41
std::vector< TrackingVertex > TrackingVertexCollection
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
T const * product() const
Definition: Handle.h:74
KVFTest(const edm::ParameterSet &)
Definition: KVFTest.cc:25
edm::EventID id() const
Definition: EventBase.h:56
edm::ParameterSet kvfPSet
Definition: KVFTest.h:57
tuple cout
Definition: gather_cfg.py:121
~KVFTest()
Definition: KVFTest.cc:37