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: " << a.get<2>() << std::endl;
81  if (recoBeamSpotHandle.isValid()) {
82  SingleTrackVertexConstraint::BTFtuple b = stvc.constrain(t_tks[i], *recoBeamSpotHandle);
83  std::cout << "Chi2: " << b.get<2>() << 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 
T getParameter(std::string const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
void beginJob() override
KVFTrackUpdate(const edm::ParameterSet &)
void endJob() override
void analyze(const edm::Event &, const edm::EventSetup &) override
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
bool isValid() const
Definition: HandleBase.h:70
double b
Definition: hdecay.h:118
BTFtuple constrain(const reco::TransientTrack &track, const GlobalPoint &priorPos, const GlobalError &priorError) const
edm::EventID id() const
Definition: EventBase.h:59
fixed size matrix
HLT enums.
double a
Definition: hdecay.h:119
T get() const
Definition: EventSetup.h:73
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
boost::tuple< bool, reco::TransientTrack, float > BTFtuple
~KVFTrackUpdate() override