CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
KVFTrackUpdate Class Reference

#include <KVFTrackUpdate.h>

Inheritance diagram for KVFTrackUpdate:
edm::one::EDAnalyzer<> edm::one::EDAnalyzerBase edm::EDConsumerBase

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void beginJob () override
 
void endJob () override
 
 KVFTrackUpdate (const edm::ParameterSet &)
 
 ~KVFTrackUpdate () override
 
- Public Member Functions inherited from edm::one::EDAnalyzer<>
 EDAnalyzer ()=default
 
 EDAnalyzer (const EDAnalyzer &)=delete
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
const EDAnalyzeroperator= (const EDAnalyzer &)=delete
 
bool wantsGlobalLuminosityBlocks () const noexcept final
 
bool wantsGlobalRuns () const noexcept final
 
bool wantsInputProcessBlocks () const noexcept final
 
bool wantsProcessBlocks () const noexcept final
 
- Public Member Functions inherited from edm::one::EDAnalyzerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const noexcept
 
bool wantsStreamRuns () const noexcept
 
 ~EDAnalyzerBase () 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
 
ESResolverIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
std::vector< ESResolverIndex > const & esGetTokenIndicesVector (edm::Transition iTrans) const
 
std::vector< ESRecordIndex > const & esGetTokenRecordIndicesVector (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::array< std::vector< ModuleDescription const *> *, NumBranchTypes > &modulesAll, std::vector< ModuleProcessName > &modulesInPreviousProcesses, 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
 
void selectInputProcessBlocks (ProductRegistry const &productRegistry, ProcessBlockHelperBase const &processBlockHelperBase)
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProductResolverIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Attributes

const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecordestoken_TTB
 
edm::EDGetTokenT< reco::BeamSpottoken_beamSpot
 
edm::EDGetTokenT< reco::TrackCollectiontoken_tracks
 

Additional Inherited Members

- Public Types inherited from edm::one::EDAnalyzerBase
typedef EDAnalyzerBase ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::one::EDAnalyzerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
template<BranchType B = InEvent>
EDConsumerBaseAdaptor< Bconsumes (edm::InputTag tag) noexcept
 
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 ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<Transition Tr = Transition::Event>
constexpr auto esConsumes ()
 
template<Transition Tr = Transition::Event>
auto esConsumes (ESInputTag tag)
 
template<Transition Tr = Transition::Event>
ESGetTokenGeneric esConsumes (eventsetup::EventSetupRecordKey const &iRecord, eventsetup::DataKey const &iKey)
 Used with EventSetupRecord::doGet. More...
 
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)
 
void resetItemsToGetFrom (BranchType iType)
 

Detailed Description

This is a very simple test analyzer to test the update of a track with a vertex constraint with the Kalman filter.

Definition at line 22 of file KVFTrackUpdate.h.

Constructor & Destructor Documentation

◆ KVFTrackUpdate()

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

Definition at line 21 of file KVFTrackUpdate.cc.

References edm::ParameterSet::getParameter(), token_beamSpot, and token_tracks.

22  : estoken_TTB(esConsumes(edm::ESInputTag("", "TransientTrackBuilder"))) {
23  token_tracks = consumes<TrackCollection>(iConfig.getParameter<InputTag>("TrackLabel"));
24  token_beamSpot = consumes<BeamSpot>(iConfig.getParameter<InputTag>("beamSpotLabel"));
25 }
edm::EDGetTokenT< reco::TrackCollection > token_tracks
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > estoken_TTB
edm::EDGetTokenT< reco::BeamSpot > token_beamSpot

◆ ~KVFTrackUpdate()

KVFTrackUpdate::~KVFTrackUpdate ( )
override

Definition at line 27 of file KVFTrackUpdate.cc.

27 {}

Member Function Documentation

◆ analyze()

void KVFTrackUpdate::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
overridevirtual

Implements edm::one::EDAnalyzerBase.

Definition at line 37 of file KVFTrackUpdate.cc.

References a, b, SingleTrackVertexConstraint::constrain(), MillePedeFileConverter_cfg::e, submitPVResolutionJobs::err, estoken_TTB, cppFunctionSkipper::exception, edm::EventSetup::getData(), mps_fire::i, iEvent, edm::HandleBase::isValid(), token_beamSpot, and token_tracks.

37  {
38  try {
39  edm::LogInfo("RecoVertex/KVFTrackUpdate") << "Reconstructing event number: " << iEvent.id() << "\n";
40 
41  // get RECO tracks from the event
42  // `tks` can be used as a ptr to a reco::TrackCollection
44  iEvent.getByToken(token_tracks, tks);
45 
46  edm::LogInfo("RecoVertex/KVFTrackUpdate") << "Found: " << (*tks).size() << " reconstructed tracks"
47  << "\n";
48  edm::LogPrint("RecoVertex/KVFTrackUpdate") << "got " << (*tks).size() << " tracks " << std::endl;
49 
50  // Transform Track to TransientTrack
51 
52  //get the builder:
53  const auto& theB = &iSetup.getData(estoken_TTB);
54  //do the conversion:
55  std::vector<TransientTrack> t_tks = theB->build(tks);
56 
57  edm::LogInfo("RecoVertex/KVFTrackUpdate") << "Found: " << t_tks.size() << " reconstructed tracks"
58  << "\n";
59 
60  GlobalPoint glbPos(0., 0., 0.);
61 
63  mat[0][0] = (20.e-04) * (20.e-04);
64  mat[1][1] = (20.e-04) * (20.e-04);
65  mat[2][2] = (5.3) * (5.3);
66  GlobalError glbErrPos(mat);
67 
68  reco::BeamSpot vertexBeamSpot;
69  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
70  iEvent.getByToken(token_beamSpot, recoBeamSpotHandle);
71 
73  for (unsigned int i = 0; i < t_tks.size(); i++) {
74  SingleTrackVertexConstraint::BTFtuple a = stvc.constrain(t_tks[i], glbPos, glbErrPos);
75  edm::LogPrint("RecoVertex/KVFTrackUpdate") << "Chi2: " << std::get<2>(a) << std::endl;
76  if (recoBeamSpotHandle.isValid()) {
77  SingleTrackVertexConstraint::BTFtuple b = stvc.constrain(t_tks[i], *recoBeamSpotHandle);
78  edm::LogPrint("RecoVertex/KVFTrackUpdate") << "Chi2: " << std::get<2>(b) << std::endl;
79  }
80  }
81  }
82 
83  catch (std::exception& err) {
84  edm::LogInfo("RecoVertex/KVFTrackUpdate") << "Exception during event number: " << iEvent.id() << "\n"
85  << err.what() << "\n";
86  }
87 }
edm::EDGetTokenT< reco::TrackCollection > token_tracks
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
std::tuple< bool, reco::TransientTrack, float > BTFtuple
BTFtuple constrain(const reco::TransientTrack &track, const GlobalPoint &priorPos, const GlobalError &priorError) const
int iEvent
Definition: GenABIO.cc:224
Log< level::Warning, true > LogPrint
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > estoken_TTB
Log< level::Info, false > LogInfo
double b
Definition: hdecay.h:120
bool isValid() const
Definition: HandleBase.h:70
double a
Definition: hdecay.h:121
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
edm::EDGetTokenT< reco::BeamSpot > token_beamSpot

◆ beginJob()

void KVFTrackUpdate::beginJob ( void  )
overridevirtual

Reimplemented from edm::one::EDAnalyzerBase.

Definition at line 29 of file KVFTrackUpdate.cc.

29 {}

◆ endJob()

void KVFTrackUpdate::endJob ( void  )
overridevirtual

Reimplemented from edm::one::EDAnalyzerBase.

Definition at line 31 of file KVFTrackUpdate.cc.

31 {}

Member Data Documentation

◆ estoken_TTB

const edm::ESGetToken<TransientTrackBuilder, TransientTrackRecord> KVFTrackUpdate::estoken_TTB
private

Definition at line 33 of file KVFTrackUpdate.h.

Referenced by analyze().

◆ token_beamSpot

edm::EDGetTokenT<reco::BeamSpot> KVFTrackUpdate::token_beamSpot
private

Definition at line 35 of file KVFTrackUpdate.h.

Referenced by analyze(), and KVFTrackUpdate().

◆ token_tracks

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

Definition at line 34 of file KVFTrackUpdate.h.

Referenced by analyze(), and KVFTrackUpdate().