CMS 3D CMS Logo

KVFTrackUpdate.cc
Go to the documentation of this file.
2 
10 
19 
20 #include <iostream>
21 
22 using namespace reco;
23 using namespace edm;
24 using namespace std;
25 
27  token_tracks = consumes<TrackCollection>(iConfig.getParameter<InputTag>("TrackLabel"));
28  token_beamSpot = consumes<BeamSpot>(iConfig.getParameter<InputTag>("beamSpotLabel"));
29 }
30 
32 
34 
36 
37 //
38 // member functions
39 //
40 
42  try {
43  edm::LogInfo("RecoVertex/KVFTrackUpdate") << "Reconstructing event number: " << iEvent.id() << "\n";
44 
45  // get RECO tracks from the event
46  // `tks` can be used as a ptr to a reco::TrackCollection
48  iEvent.getByToken(token_tracks, tks);
49 
50  edm::LogInfo("RecoVertex/KVFTrackUpdate") << "Found: " << (*tks).size() << " reconstructed tracks"
51  << "\n";
52  std::cout << "got " << (*tks).size() << " tracks " << std::endl;
53 
54  // Transform Track to TransientTrack
55 
56  //get the builder:
58  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder", theB);
59  //do the conversion:
60  std::vector<TransientTrack> t_tks = (*theB).build(tks);
61 
62  edm::LogInfo("RecoVertex/KVFTrackUpdate") << "Found: " << t_tks.size() << " reconstructed tracks"
63  << "\n";
64 
65  GlobalPoint glbPos(0., 0., 0.);
66 
68  mat[0][0] = (20.e-04) * (20.e-04);
69  mat[1][1] = (20.e-04) * (20.e-04);
70  mat[2][2] = (5.3) * (5.3);
71  GlobalError glbErrPos(mat);
72 
73  reco::BeamSpot vertexBeamSpot;
74  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
75  iEvent.getByToken(token_beamSpot, recoBeamSpotHandle);
76 
78  for (unsigned int i = 0; i < t_tks.size(); i++) {
79  SingleTrackVertexConstraint::BTFtuple a = stvc.constrain(t_tks[i], glbPos, glbErrPos);
80  std::cout << "Chi2: " << std::get<2>(a) << std::endl;
81  if (recoBeamSpotHandle.isValid()) {
82  SingleTrackVertexConstraint::BTFtuple b = stvc.constrain(t_tks[i], *recoBeamSpotHandle);
83  std::cout << "Chi2: " << std::get<2>(b) << std::endl;
84  }
85  }
86  }
87 
88  catch (std::exception& err) {
89  edm::LogInfo("RecoVertex/KVFTrackUpdate") << "Exception during event number: " << iEvent.id() << "\n"
90  << err.what() << "\n";
91  }
92 }
93 
Handle.h
AlgebraicSymMatrix33
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
Definition: AlgebraicROOTObjects.h:21
mps_fire.i
i
Definition: mps_fire.py:428
MessageLogger.h
KalmanVertexFitter.h
ESHandle.h
edm
HLT enums.
Definition: AlignableModifier.h:19
gather_cfg.cout
cout
Definition: gather_cfg.py:144
KVFTrackUpdate::beginJob
void beginJob() override
Definition: KVFTrackUpdate.cc:33
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
edm::LogInfo
Log< level::Info, false > LogInfo
Definition: MessageLogger.h:125
TransientTrack.h
edm::Handle< reco::TrackCollection >
MakerMacros.h
SingleTrackVertexConstraint.h
Track.h
edm::EventSetup::get
T get() const
Definition: EventSetup.h:80
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
SingleTrackVertexConstraint::constrain
BTFtuple constrain(const reco::TransientTrack &track, const GlobalPoint &priorPos, const GlobalError &priorError) const
Definition: SingleTrackVertexConstraint.cc:21
reco::BeamSpot
Definition: BeamSpot.h:21
TransientTrackRecord
Definition: TransientTrackRecord.h:11
edm::ESHandle< TransientTrackBuilder >
SingleTrackVertexConstraint::BTFtuple
std::tuple< bool, reco::TransientTrack, float > BTFtuple
Definition: SingleTrackVertexConstraint.h:23
Point3DBase< float, GlobalTag >
b
double b
Definition: hdecay.h:118
Vertex.h
cppFunctionSkipper.exception
exception
Definition: cppFunctionSkipper.py:10
KVFTrackUpdate.h
KVFTrackUpdate::endJob
void endJob() override
Definition: KVFTrackUpdate.cc:35
TransientTrackBuilder.h
edm::ParameterSet
Definition: ParameterSet.h:47
a
double a
Definition: hdecay.h:119
KVFTrackUpdate
Definition: KVFTrackUpdate.h:20
iEvent
int iEvent
Definition: GenABIO.cc:224
GlobalErrorBase< double, ErrorMatrixTag >
edm::EventSetup
Definition: EventSetup.h:57
TrackAssociatorRecord.h
TransientTrackRecord.h
submitPVResolutionJobs.err
err
Definition: submitPVResolutionJobs.py:85
InputTag.h
reco::get
T get(const Candidate &c)
Definition: component.h:60
VertexFwd.h
TransientVertex.h
std
Definition: JetResolutionObject.h:76
KVFTrackUpdate::~KVFTrackUpdate
~KVFTrackUpdate() override
Definition: KVFTrackUpdate.cc:31
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::HandleBase::isValid
bool isValid() const
Definition: HandleBase.h:70
KVFTrackUpdate::KVFTrackUpdate
KVFTrackUpdate(const edm::ParameterSet &)
Definition: KVFTrackUpdate.cc:26
edm::Event
Definition: Event.h:73
edm::InputTag
Definition: InputTag.h:15
KVFTrackUpdate::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: KVFTrackUpdate.cc:41
SingleTrackVertexConstraint
Definition: SingleTrackVertexConstraint.h:20
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37