CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
KVFTest Class Reference

#include <RecoVertex/KVFTest/src/KVFTest.cc>

Inheritance diagram for KVFTest:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
 
virtual void beginJob ()
 
virtual void endJob ()
 
 KVFTest (const edm::ParameterSet &)
 
 ~KVFTest ()
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

TrackingVertex getSimVertex (const edm::Event &iEvent) const
 

Private Attributes

TrackAssociatorByChi2associatorForParamAtPca
 
edm::ParameterSet kvfPSet
 
std::string outputFile_
 
TFile * rootFile_
 
edm::ParameterSet theConfig
 
edm::EDGetTokenT
< reco::TrackCollection
token_tracks
 
edm::EDGetTokenT
< TrackingParticleCollection
token_TrackTruth
 
edm::EDGetTokenT
< TrackingVertexCollection
token_VertexTruth
 
SimpleVertexTreetree
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

This is a very simple test analyzer mean to test the KalmanVertexFitter

Description: steers tracker primary vertex reconstruction and storage

Implementation: <Notes on="" implementation>="">

Definition at line 37 of file KVFTest.h.

Constructor & Destructor Documentation

KVFTest::KVFTest ( const edm::ParameterSet iConfig)
explicit

Definition at line 24 of file KVFTest.cc.

References edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), kvfPSet, outputFile_, rootFile_, AlCaHLTBitMon_QueryRunRegistry::string, token_tracks, token_TrackTruth, and token_VertexTruth.

25  : theConfig(iConfig), associatorForParamAtPca(0), tree(0)
26 {
27  token_tracks = consumes<TrackCollection>(iConfig.getParameter<string>("TrackLabel"));
28  outputFile_ = iConfig.getUntrackedParameter<std::string>("outputFile");
29  kvfPSet = iConfig.getParameter<edm::ParameterSet>("KVFParameters");
30  rootFile_ = TFile::Open(outputFile_.c_str(),"RECREATE");
31  edm::LogInfo("RecoVertex/KVFTest")
32  << "Initializing KVF TEST analyser - Output file: " << outputFile_ <<"\n";
33 
34  token_TrackTruth = consumes<TrackingParticleCollection>(edm::InputTag("trackingtruth", "TrackTruth"));
35  token_VertexTruth = consumes<TrackingVertexCollection>(edm::InputTag("trackingtruth", "VertexTruth"));
36 
37 }
TrackAssociatorByChi2 * associatorForParamAtPca
Definition: KVFTest.h:53
T getParameter(std::string const &) const
std::string outputFile_
Definition: KVFTest.h:57
T getUntrackedParameter(std::string const &, T const &) const
SimpleVertexTree * tree
Definition: KVFTest.h:54
edm::EDGetTokenT< reco::TrackCollection > token_tracks
Definition: KVFTest.h:58
edm::ParameterSet theConfig
Definition: KVFTest.h:51
edm::EDGetTokenT< TrackingParticleCollection > token_TrackTruth
Definition: KVFTest.h:59
TFile * rootFile_
Definition: KVFTest.h:55
edm::EDGetTokenT< TrackingVertexCollection > token_VertexTruth
Definition: KVFTest.h:60
edm::ParameterSet kvfPSet
Definition: KVFTest.h:52
KVFTest::~KVFTest ( )

Definition at line 40 of file KVFTest.cc.

References rootFile_.

40  {
41  delete rootFile_;
42 }
TFile * rootFile_
Definition: KVFTest.h:55

Member Function Documentation

void KVFTest::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
virtual

Implements edm::EDAnalyzer.

Definition at line 57 of file KVFTest.cc.

References TrackAssociatorByChi2::associateRecoToSim(), associatorForParamAtPca, gather_cfg::cout, SimpleVertexTree::fill(), edm::EventSetup::get(), edm::Event::getByToken(), getSimVertex(), edm::EventBase::id(), edm::HandleBase::isValid(), TransientVertex::position(), edm::Handle< T >::product(), edm::ESHandle< class >::product(), token_tracks, token_TrackTruth, tree, and KalmanVertexFitter::vertex().

58 {
59  if ( associatorForParamAtPca==0 ) {
60  edm::ESHandle<TrackAssociatorBase> theAssociatorForParamAtPca;
61  iSetup.get<TrackAssociatorRecord>().get("TrackAssociatorByChi2",theAssociatorForParamAtPca);
62  associatorForParamAtPca = (TrackAssociatorByChi2 *) theAssociatorForParamAtPca.product();
63 
64  tree = new SimpleVertexTree("VertexFitter", associatorForParamAtPca);
65  }
66 
67 
68 
69  edm::LogInfo("RecoVertex/KVFTest")
70  << "Reconstructing event number: " << iEvent.id() << "\n";
71 
72  // get RECO tracks from the event
73  // `tks` can be used as a ptr to a reco::TrackCollection
75  iEvent.getByToken(token_tracks, tks);
76  if (!tks.isValid()) {
77  edm::LogInfo("RecoVertex/KVFTest")
78  << "Exception during event number: " << iEvent.id()
79  << "\n";
80  } else {
81  edm::LogInfo("RecoVertex/KVFTest")
82  << "Found: " << (*tks).size() << " reconstructed tracks" << "\n";
83  std::cout << "got " << (*tks).size() << " tracks " << std::endl;
84 
85  // Transform Track to TransientTrack
86 
87  //get the builder:
89  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
90  //do the conversion:
91  std::vector<TransientTrack> t_tks = (*theB).build(tks);
92 
93  edm::LogInfo("RecoVertex/KVFTest")
94  << "Found: " << t_tks.size() << " reconstructed tracks" << "\n";
95 
96  // Call the KalmanVertexFitter if more than 1 track
97  if (t_tks.size() > 1) {
98  // KalmanVertexFitter kvf(kvfPSet);
99  KalmanVertexFitter kvf(true);
100  TransientVertex tv = kvf.vertex(t_tks);
101 
102  std::cout << "Position: " << Vertex::Point(tv.position()) << "\n";
103 
104  // For the analysis: compare to your SimVertex
105  TrackingVertex sv = getSimVertex(iEvent);
107  iEvent.getByToken(token_TrackTruth, TPCollectionH);
108  if (!TPCollectionH.isValid()) {
109  edm::LogInfo("RecoVertex/KVFTest")
110  << "Exception during event number: " << iEvent.id()
111  << "\n";
112  } else {
113  const TrackingParticleCollection tPC = *(TPCollectionH.product());
115  TPCollectionH,
116  &iEvent);
117  tree->fill(tv, &sv, &recSimColl);
118  }
119  }
120  }
121 }
TrackAssociatorByChi2 * associatorForParamAtPca
Definition: KVFTest.h:53
std::vector< TrackingParticle > TrackingParticleCollection
TrackingVertex getSimVertex(const edm::Event &iEvent) const
Definition: KVFTest.cc:125
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
SimpleVertexTree * tree
Definition: KVFTest.h:54
edm::EDGetTokenT< reco::TrackCollection > token_tracks
Definition: KVFTest.h:58
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
edm::EDGetTokenT< TrackingParticleCollection > token_TrackTruth
Definition: KVFTest.h:59
GlobalPoint position() const
bool isValid() const
Definition: HandleBase.h:76
double associateRecoToSim(reco::TrackCollection::const_iterator, TrackingParticleCollection::const_iterator, const reco::BeamSpot &) const
compare reco::TrackCollection and TrackingParticleCollection iterators: returns the chi2 ...
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
T const * product() const
Definition: Handle.h:81
edm::EventID id() const
Definition: EventBase.h:56
tuple cout
Definition: gather_cfg.py:121
void KVFTest::beginJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 44 of file KVFTest.cc.

44  {
45 }
void KVFTest::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 48 of file KVFTest.cc.

References tree.

48  {
49  delete tree;
50 }
SimpleVertexTree * tree
Definition: KVFTest.h:54
TrackingVertex KVFTest::getSimVertex ( const edm::Event iEvent) const
private

Definition at line 125 of file KVFTest.cc.

References edm::Event::getByToken(), edm::Handle< T >::product(), and token_VertexTruth.

Referenced by analyze().

126 {
127  // get the simulated vertices
129  iEvent.getByToken(token_VertexTruth,TVCollectionH);
130  const TrackingVertexCollection tPC = *(TVCollectionH.product());
131 
132 // Handle<edm::SimVertexContainer> simVtcs;
133 // iEvent.getByLabel("g4SimHits", simVtcs);
134 // std::cout << "SimVertex " << simVtcs->size() << std::endl;
135 // for(edm::SimVertexContainer::const_iterator v=simVtcs->begin();
136 // v!=simVtcs->end(); ++v){
137 // std::cout << "simvtx "
138 // << v->position().x() << " "
139 // << v->position().y() << " "
140 // << v->position().z() << " "
141 // << v->parentIndex() << " "
142 // << v->noParent() << " "
143 // << std::endl;
144 // }
145  return *(tPC.begin());
146 }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
edm::EDGetTokenT< TrackingVertexCollection > token_VertexTruth
Definition: KVFTest.h:60
std::vector< TrackingVertex > TrackingVertexCollection
T const * product() const
Definition: Handle.h:81

Member Data Documentation

TrackAssociatorByChi2* KVFTest::associatorForParamAtPca
private

Definition at line 53 of file KVFTest.h.

Referenced by analyze().

edm::ParameterSet KVFTest::kvfPSet
private

Definition at line 52 of file KVFTest.h.

Referenced by KVFTest().

std::string KVFTest::outputFile_
private

Definition at line 57 of file KVFTest.h.

Referenced by KVFTest().

TFile* KVFTest::rootFile_
private

Definition at line 55 of file KVFTest.h.

Referenced by KVFTest(), and ~KVFTest().

edm::ParameterSet KVFTest::theConfig
private

Definition at line 51 of file KVFTest.h.

edm::EDGetTokenT<reco::TrackCollection> KVFTest::token_tracks
private

Definition at line 58 of file KVFTest.h.

Referenced by analyze(), and KVFTest().

edm::EDGetTokenT<TrackingParticleCollection> KVFTest::token_TrackTruth
private

Definition at line 59 of file KVFTest.h.

Referenced by analyze(), and KVFTest().

edm::EDGetTokenT<TrackingVertexCollection> KVFTest::token_VertexTruth
private

Definition at line 60 of file KVFTest.h.

Referenced by getSimVertex(), and KVFTest().

SimpleVertexTree* KVFTest::tree
private