CMS 3D CMS Logo

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

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void beginJob () override
 
void endJob () override
 
 KVFTest (const edm::ParameterSet &)
 
 ~KVFTest () override
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
SerialTaskQueueglobalLuminosityBlocksQueue ()
 
SerialTaskQueueglobalRunsQueue ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
 ~EDAnalyzer () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

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

Private Attributes

edm::ParameterSet kvfPSet
 
std::string outputFile_
 
TFile * rootFile_
 
edm::ParameterSet theConfig
 
edm::EDGetTokenT< reco::TrackToTrackingParticleAssociatortoken_associatorForParamAtPca
 
edm::EDGetTokenT< reco::TrackCollectiontoken_tracks
 
edm::EDGetTokenT< TrackingParticleCollectiontoken_TrackTruth
 
edm::EDGetTokenT< TrackingVertexCollectiontoken_VertexTruth
 
std::unique_ptr< SimpleVertexTreetree
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
static bool wantsGlobalLuminosityBlocks ()
 
static bool wantsGlobalRuns ()
 
static bool wantsStreamLuminosityBlocks ()
 
static bool wantsStreamRuns ()
 
- 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 ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
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 36 of file KVFTest.h.

Constructor & Destructor Documentation

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

Definition at line 25 of file KVFTest.cc.

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

25  : 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 }
T getParameter(std::string const &) const
std::string outputFile_
Definition: KVFTest.h:54
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< reco::TrackCollection > token_tracks
Definition: KVFTest.h:55
edm::ParameterSet theConfig
Definition: KVFTest.h:49
edm::EDGetTokenT< reco::TrackToTrackingParticleAssociator > token_associatorForParamAtPca
Definition: KVFTest.h:58
edm::EDGetTokenT< TrackingParticleCollection > token_TrackTruth
Definition: KVFTest.h:56
TFile * rootFile_
Definition: KVFTest.h:52
edm::EDGetTokenT< TrackingVertexCollection > token_VertexTruth
Definition: KVFTest.h:57
edm::ParameterSet kvfPSet
Definition: KVFTest.h:50
KVFTest::~KVFTest ( )
override

Definition at line 38 of file KVFTest.cc.

References rootFile_.

38 { delete rootFile_; }
TFile * rootFile_
Definition: KVFTest.h:52

Member Function Documentation

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

Definition at line 48 of file KVFTest.cc.

References reco::TrackToTrackingParticleAssociator::associateRecoToSim(), gather_cfg::cout, edm::EventSetup::get(), edm::Event::getByToken(), getSimVertex(), edm::EventBase::id(), edm::HandleBase::isValid(), TransientVertex::position(), edm::Handle< T >::product(), edm::ESHandle< T >::product(), pfDeepBoostedJetPreprocessParams_cfi::sv, token_associatorForParamAtPca, token_tracks, token_TrackTruth, and KalmanVertexFitter::vertex().

48  {
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
91  TrackingVertex sv = getSimVertex(iEvent);
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 }
std::vector< TrackingParticle > TrackingParticleCollection
TrackingVertex getSimVertex(const edm::Event &iEvent) const
Definition: KVFTest.cc:107
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
edm::EDGetTokenT< reco::TrackCollection > token_tracks
Definition: KVFTest.h:55
edm::EDGetTokenT< reco::TrackToTrackingParticleAssociator > token_associatorForParamAtPca
Definition: KVFTest.h:58
edm::EDGetTokenT< TrackingParticleCollection > token_TrackTruth
Definition: KVFTest.h:56
GlobalPoint position() const
math::XYZPoint Point
bool isValid() const
Definition: HandleBase.h:70
reco::RecoToSimCollection associateRecoToSim(const edm::Handle< edm::View< reco::Track >> &tCH, const edm::Handle< TrackingParticleCollection > &tPCH) const
T const * product() const
Definition: Handle.h:69
edm::EventID id() const
Definition: EventBase.h:59
T get() const
Definition: EventSetup.h:73
Definition: tree.py:1
T const * product() const
Definition: ESHandle.h:86
void KVFTest::beginJob ( void  )
overridevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 40 of file KVFTest.cc.

40 {}
void KVFTest::endJob ( void  )
overridevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 42 of file KVFTest.cc.

42 {}
TrackingVertex KVFTest::getSimVertex ( const edm::Event iEvent) const
private

Definition at line 107 of file KVFTest.cc.

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

Referenced by analyze().

107  {
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 }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
T const * product() const
Definition: Handle.h:69
edm::EDGetTokenT< TrackingVertexCollection > token_VertexTruth
Definition: KVFTest.h:57
std::vector< TrackingVertex > TrackingVertexCollection

Member Data Documentation

edm::ParameterSet KVFTest::kvfPSet
private

Definition at line 50 of file KVFTest.h.

Referenced by KVFTest().

std::string KVFTest::outputFile_
private

Definition at line 54 of file KVFTest.h.

Referenced by KVFTest().

TFile* KVFTest::rootFile_
private

Definition at line 52 of file KVFTest.h.

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

edm::ParameterSet KVFTest::theConfig
private

Definition at line 49 of file KVFTest.h.

edm::EDGetTokenT<reco::TrackToTrackingParticleAssociator> KVFTest::token_associatorForParamAtPca
private

Definition at line 58 of file KVFTest.h.

Referenced by analyze(), and KVFTest().

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

Definition at line 55 of file KVFTest.h.

Referenced by analyze(), and KVFTest().

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

Definition at line 56 of file KVFTest.h.

Referenced by analyze(), and KVFTest().

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

Definition at line 57 of file KVFTest.h.

Referenced by getSimVertex(), and KVFTest().

std::unique_ptr<SimpleVertexTree> KVFTest::tree
private