CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackAlgoCompareUtil.cc
Go to the documentation of this file.
2 
3 using namespace std;
4 using namespace edm;
5 
6 
7 // constructors and destructor
9  trackLabel_algoA(consumes<View<reco::Track>>(iConfig.getParameter<edm::InputTag>("trackLabel_algoA"))),
10  trackLabel_algoB(consumes<View<reco::Track>>(iConfig.getParameter<edm::InputTag>("trackLabel_algoB"))),
11  trackingParticleLabel_fakes(consumes<TrackingParticleCollection>(iConfig.getParameter<edm::InputTag>("trackingParticleLabel_fakes"))),
12  trackingParticleLabel_effic(consumes<TrackingParticleCollection>(iConfig.getParameter<edm::InputTag>("trackingParticleLabel_effic"))),
13  beamSpotLabel(consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpotLabel"))),
14  UseAssociators(iConfig.getParameter< bool >("UseAssociators")),
15  UseVertex(iConfig.getParameter< bool >("UseVertex"))
16 {
17  //now do what ever other initialization is needed
18  if(UseVertex) {
19  vertexLabel_algoA = consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertexLabel_algoA"));
20  vertexLabel_algoB = consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertexLabel_algoB"));
21  }
22 
23  if(UseAssociators) {
24  assocLabel_algoA = consumes<reco::TrackToTrackingParticleAssociator>(iConfig.getUntrackedParameter<std::string>("assocLabel_algoA", "trackAssociatorByHits"));
25  assocLabel_algoB = consumes<reco::TrackToTrackingParticleAssociator>(iConfig.getUntrackedParameter<std::string>("assocLabel_algoB", "trackAssociatorByHits"));
26  }
27  else {
28  edm::InputTag algoA = iConfig.getParameter< edm::InputTag >("associatormap_algoA");
29  edm::InputTag algoB = iConfig.getParameter< edm::InputTag >("associatormap_algoB");
30 
31  associatormap_algoA_recoToSim = consumes<reco::RecoToSimCollection>(algoA);
32  associatormap_algoB_recoToSim = consumes<reco::RecoToSimCollection>(algoB);
33  associatormap_algoA_simToReco = consumes<reco::SimToRecoCollection>(algoA);
34  associatormap_algoB_simToReco = consumes<reco::SimToRecoCollection>(algoB);
35  }
36 
37  produces<RecoTracktoTPCollection>("AlgoA");
38  produces<RecoTracktoTPCollection>("AlgoB");
39  produces<TPtoRecoTrackCollection>("TP");
40 }
41 
42 
44 {
45 }
46 
47 
48 // ------------ method called to produce the data ------------
49 void
51 {
52  // create output collection instance
53  std::auto_ptr<RecoTracktoTPCollection> outputAlgoA(new RecoTracktoTPCollection());
54  std::auto_ptr<RecoTracktoTPCollection> outputAlgoB(new RecoTracktoTPCollection());
55  std::auto_ptr<TPtoRecoTrackCollection> outputTP(new TPtoRecoTrackCollection());
56 
57  // Get Inputs
58  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
59  iEvent.getByToken(beamSpotLabel, recoBeamSpotHandle);
60  reco::BeamSpot beamSpot = *recoBeamSpotHandle;
61 
62  edm::Handle<View<reco::Track> > trackCollAlgoA;
63  iEvent.getByToken(trackLabel_algoA, trackCollAlgoA);
64 
65  edm::Handle< View<reco::Track> > trackCollAlgoB;
66  iEvent.getByToken(trackLabel_algoB, trackCollAlgoB);
67 
68  edm::Handle<TrackingParticleCollection> trackingParticleCollFakes;
69  iEvent.getByToken(trackingParticleLabel_fakes, trackingParticleCollFakes);
70 
71  edm::Handle<TrackingParticleCollection> trackingParticleCollEffic;
72  iEvent.getByToken(trackingParticleLabel_effic, trackingParticleCollEffic);
73 
76  if(UseVertex)
77  {
78  iEvent.getByToken(vertexLabel_algoA, vertexCollAlgoA);
79  iEvent.getByToken(vertexLabel_algoB, vertexCollAlgoB);
80  }
81 
82  // call the associator functions:
83  reco::RecoToSimCollection recSimColl_AlgoA;
84  reco::RecoToSimCollection recSimColl_AlgoB;
85 
86  reco::SimToRecoCollection simRecColl_AlgoA;
87  reco::SimToRecoCollection simRecColl_AlgoB;
88 
89  if(UseAssociators)
90  {
92  iEvent.getByToken(assocLabel_algoA, theAssociator_algoA);
93 
95  iEvent.getByToken(assocLabel_algoB, theAssociator_algoB);
96 
97  recSimColl_AlgoA = theAssociator_algoA->associateRecoToSim(trackCollAlgoA, trackingParticleCollFakes);
98  recSimColl_AlgoB = theAssociator_algoB->associateRecoToSim(trackCollAlgoB, trackingParticleCollFakes);
99 
100  simRecColl_AlgoA = theAssociator_algoA->associateSimToReco(trackCollAlgoA, trackingParticleCollEffic);
101  simRecColl_AlgoB = theAssociator_algoB->associateSimToReco(trackCollAlgoB, trackingParticleCollEffic);
102  }
103  else
104  {
105  Handle<reco::RecoToSimCollection > recotosimCollectionH_AlgoA;
106  iEvent.getByToken(associatormap_algoA_recoToSim,recotosimCollectionH_AlgoA);
107  recSimColl_AlgoA = *(recotosimCollectionH_AlgoA.product());
108 
109  Handle<reco::RecoToSimCollection > recotosimCollectionH_AlgoB;
110  iEvent.getByToken(associatormap_algoB_recoToSim,recotosimCollectionH_AlgoB);
111  recSimColl_AlgoB = *(recotosimCollectionH_AlgoB.product());
112 
113  Handle<reco::SimToRecoCollection > simtorecoCollectionH_AlgoA;
114  iEvent.getByToken(associatormap_algoA_simToReco, simtorecoCollectionH_AlgoA);
115  simRecColl_AlgoA = *(simtorecoCollectionH_AlgoA.product());
116 
117  Handle<reco::SimToRecoCollection > simtorecoCollectionH_AlgoB;
118  iEvent.getByToken(associatormap_algoB_simToReco, simtorecoCollectionH_AlgoB);
119  simRecColl_AlgoB = *(simtorecoCollectionH_AlgoB.product());
120  }
121 
122  // define the vector of references to trackingParticleColl associated with a given reco::Track
123  std::vector<std::pair<TrackingParticleRef, double> > associatedTrackingParticles;
124 
125  // define the vector of references to trackColl associated with a given TrackingParticle
126  std::vector<std::pair<reco::TrackBaseRef, double> > associatedRecoTracks;
127 
128  // Get the magnetic field data from the event (used to calculate the point of closest TrackingParticle)
129  edm::ESHandle<MagneticField> theMagneticField;
130  iSetup.get<IdealMagneticFieldRecord>().get(theMagneticField);
131  const MagneticField *magneticField = theMagneticField.product();
132 
133  // fill collection algoA
134  for(View<reco::Track>::size_type i = 0; i < trackCollAlgoA->size(); ++i)
135  {
136  // get recoTrack algo A
137  reco::TrackBaseRef recoTrack(trackCollAlgoA, i);
138  RecoTracktoTP recoTracktoTP;
139  recoTracktoTP.SetRecoTrack(recoTrack);
140  recoTracktoTP.SetBeamSpot(beamSpot.position());
141 
142  // get the associated trackingParticle
143  if(recSimColl_AlgoA.find(recoTrack) != recSimColl_AlgoA.end())
144  {
145  associatedTrackingParticles = recSimColl_AlgoA[recoTrack];
146  recoTracktoTP.SetTrackingParticle( associatedTrackingParticles.begin()->first );
147  recoTracktoTP.SetShared( associatedTrackingParticles.begin()->second );
148  SetTrackingParticleD0Dz(associatedTrackingParticles.begin()->first, beamSpot, magneticField, recoTracktoTP);
149  }
150  else
151  {
152  recoTracktoTP.SetTrackingParticle(TrackingParticleRef());
153  recoTracktoTP.SetShared(-1.0);
154  }
155 
156  // get the reco primary vertex info
157  if(UseVertex && vertexCollAlgoA->size())
158  {
159  recoTracktoTP.SetRecoVertex( reco::VertexRef(vertexCollAlgoA, 0) );
160  }
161  else
162  {
163  recoTracktoTP.SetRecoVertex( reco::VertexRef() );
164  }
165 
166  outputAlgoA->push_back(recoTracktoTP);
167  }
168 
169 
170  // fill collection algoB
171  for(reco::TrackCollection::size_type i = 0; i < trackCollAlgoB->size(); ++i)
172  {
173  // get recoTrack algo B
174  reco::TrackBaseRef recoTrack(trackCollAlgoB, i);
175  RecoTracktoTP recoTracktoTP;
176  recoTracktoTP.SetRecoTrack(recoTrack);
177  recoTracktoTP.SetBeamSpot(beamSpot.position());
178 
179  // get the associated trackingParticle
180  if(recSimColl_AlgoB.find(recoTrack) != recSimColl_AlgoB.end())
181  {
182  associatedTrackingParticles = recSimColl_AlgoB[recoTrack];
183  recoTracktoTP.SetTrackingParticle( associatedTrackingParticles.begin()->first );
184  recoTracktoTP.SetShared( associatedTrackingParticles.begin()->second );
185  SetTrackingParticleD0Dz(associatedTrackingParticles.begin()->first, beamSpot, magneticField, recoTracktoTP);
186  }
187  else
188  {
189  recoTracktoTP.SetTrackingParticle(TrackingParticleRef());
190  recoTracktoTP.SetShared(-1.0);
191  }
192 
193  // get the reco primary vertex info
194  if(UseVertex && vertexCollAlgoB->size())
195  {
196  recoTracktoTP.SetRecoVertex( reco::VertexRef(vertexCollAlgoB, 0) );
197  }
198  else
199  {
200  recoTracktoTP.SetRecoVertex( reco::VertexRef() );
201  }
202 
203  outputAlgoB->push_back(recoTracktoTP);
204  }
205 
206 
207  for(TrackingParticleCollection::size_type i = 0; i < trackingParticleCollEffic->size(); ++i)
208  {
209  // initialize the trackingParticle (sim) info
210  TrackingParticleRef tparticle(trackingParticleCollEffic, i);
211  TPtoRecoTrack tptoRecoTrack;
212  tptoRecoTrack.SetBeamSpot(beamSpot.position());
213  tptoRecoTrack.SetTrackingParticle(tparticle);
214  SetTrackingParticleD0Dz(tparticle, beamSpot, magneticField, tptoRecoTrack);
215 
216  // get the assocated recoTrack algoA
217  if(simRecColl_AlgoA.find(tparticle) != simRecColl_AlgoA.end())
218  {
219  associatedRecoTracks = simRecColl_AlgoA[tparticle];
220  tptoRecoTrack.SetRecoTrack_AlgoA(associatedRecoTracks.begin()->first );
221  tptoRecoTrack.SetShared_AlgoA(associatedRecoTracks.begin()->second );
222  }
223  else
224  {
225  tptoRecoTrack.SetRecoTrack_AlgoA(reco::TrackBaseRef());
226  tptoRecoTrack.SetShared_AlgoA(-1.0);
227  }
228 
229  // get the recoVertex algo A
230  if(UseVertex && vertexCollAlgoA->size())
231  {
232  tptoRecoTrack.SetRecoVertex_AlgoA( reco::VertexRef(vertexCollAlgoA, 0) );
233  }
234  else
235  {
236  tptoRecoTrack.SetRecoVertex_AlgoA( reco::VertexRef() );
237  }
238 
239  // get the assocated recoTrack algoB
240  if(simRecColl_AlgoB.find(tparticle) != simRecColl_AlgoB.end())
241  {
242  associatedRecoTracks = simRecColl_AlgoB[tparticle];
243  tptoRecoTrack.SetRecoTrack_AlgoB(associatedRecoTracks.begin()->first );
244  tptoRecoTrack.SetShared_AlgoB(associatedRecoTracks.begin()->second );
245  }
246  else
247  {
248  tptoRecoTrack.SetRecoTrack_AlgoB(reco::TrackBaseRef());
249  tptoRecoTrack.SetShared_AlgoB(-1.0);
250  }
251  // get the recoVertex algo B
252  if(UseVertex && vertexCollAlgoB->size())
253  {
254  tptoRecoTrack.SetRecoVertex_AlgoB( reco::VertexRef(vertexCollAlgoB, 0) );
255  }
256  else
257  {
258  tptoRecoTrack.SetRecoVertex_AlgoB( reco::VertexRef() );
259  }
260 
261  outputTP->push_back(tptoRecoTrack);
262  }
263 
264 
265  // put the collection in the event record
266  iEvent.put(outputAlgoA, "AlgoA");
267  iEvent.put(outputAlgoB, "AlgoB");
268  iEvent.put(outputTP, "TP");
269 }
270 
271 // ------------ Producer Specific Meber Fucntions ----------------------------------------
273 {
274  GlobalPoint trackingParticleVertex( tp->vertex().x(), tp->vertex().y(), tp->vertex().z() );
275  GlobalVector trackingParticleP3(tp->g4Track_begin()->momentum().x(),
276  tp->g4Track_begin()->momentum().y(),
277  tp->g4Track_begin()->momentum().z() );
278  TrackCharge trackingParticleCharge(tp->charge());
279 
280  FreeTrajectoryState ftsAtProduction( trackingParticleVertex, trackingParticleP3, trackingParticleCharge, bf );
281  TSCBLBuilderNoMaterial tscblBuilder;
282  TrajectoryStateClosestToBeamLine tsAtClosestApproach = tscblBuilder(ftsAtProduction, bs); //as in TrackProducerAlgorithm
283 
284  if(tsAtClosestApproach.isValid())
285  {
286  GlobalPoint v1 = tsAtClosestApproach.trackStateAtPCA().position();
287  GlobalVector p = tsAtClosestApproach.trackStateAtPCA().momentum();
288 
290  TPRT.SetTrackingParticlePCA(v1);
291  }
292  else
293  {
294  TPRT.SetTrackingParticleMomentumPCA(GlobalVector(-9999.0, -9999.0, -9999.0));
295  TPRT.SetTrackingParticlePCA(GlobalPoint(-9999.0, -9999.0, -9999.0));
296  }
297 }
298 
299 
301 {
302  GlobalPoint trackingParticleVertex( tp->vertex().x(), tp->vertex().y(), tp->vertex().z() );
303  GlobalVector trackingParticleP3(tp->g4Track_begin()->momentum().x(),
304  tp->g4Track_begin()->momentum().y(),
305  tp->g4Track_begin()->momentum().z() );
306  TrackCharge trackingParticleCharge(tp->charge());
307 
308  FreeTrajectoryState ftsAtProduction( trackingParticleVertex, trackingParticleP3, trackingParticleCharge, bf );
309  TSCBLBuilderNoMaterial tscblBuilder;
310  TrajectoryStateClosestToBeamLine tsAtClosestApproach = tscblBuilder(ftsAtProduction, bs); //as in TrackProducerAlgorithm
311 
312  if(tsAtClosestApproach.isValid())
313  {
314  GlobalPoint v1 = tsAtClosestApproach.trackStateAtPCA().position();
315  GlobalVector p = tsAtClosestApproach.trackStateAtPCA().momentum();
316 
318  RTTP.SetTrackingParticlePCA(v1);
319  }
320  else
321  {
322  RTTP.SetTrackingParticleMomentumPCA(GlobalVector(-9999.0, -9999.0, -9999.0));
323  RTTP.SetTrackingParticlePCA(GlobalPoint(-9999.0, -9999.0, -9999.0));
324  }
325 }
326 
T getParameter(std::string const &) const
unsigned int size_type
Definition: View.h:89
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
void SetTrackingParticle(TrackingParticleRef tp)
Definition: RecoTracktoTP.h:28
void SetTrackingParticleMomentumPCA(const GlobalVector &p)
Definition: RecoTracktoTP.h:56
std::vector< TrackingParticle > TrackingParticleCollection
const_iterator end() const
last iterator over the map (read only)
void SetTrackingParticle(TrackingParticleRef tp)
Definition: TPtoRecoTrack.h:28
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
void SetShared_AlgoA(const float &mA)
Definition: TPtoRecoTrack.h:33
void SetTrackingParticleMomentumPCA(const GlobalVector &p)
Definition: TPtoRecoTrack.h:81
edm::EDGetTokenT< reco::SimToRecoCollection > associatormap_algoB_simToReco
edm::EDGetTokenT< reco::TrackToTrackingParticleAssociator > assocLabel_algoA
const_iterator find(const key_type &k) const
find element with specified reference key
void SetShared(const float &m)
Definition: RecoTracktoTP.h:32
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
edm::EDGetTokenT< reco::VertexCollection > vertexLabel_algoA
void SetTrackingParticlePCA(const GlobalPoint &v)
Definition: RecoTracktoTP.h:57
std::vector< RecoTracktoTP > RecoTracktoTPCollection
edm::EDGetTokenT< TrackingParticleCollection > trackingParticleLabel_effic
void SetRecoVertex(reco::VertexRef vertex)
Definition: RecoTracktoTP.h:30
std::vector< TPtoRecoTrack > TPtoRecoTrackCollection
void SetRecoTrack_AlgoA(reco::TrackBaseRef track)
Definition: TPtoRecoTrack.h:30
void SetTrackingParticleD0Dz(TrackingParticleRef tp, const reco::BeamSpot &bs, const MagneticField *bf, TPtoRecoTrack &TPRT) const
void SetTrackingParticlePCA(const GlobalPoint &v)
Definition: TPtoRecoTrack.h:82
uint16_t size_type
int TrackCharge
Definition: TrackCharge.h:4
edm::EDGetTokenT< reco::RecoToSimCollection > associatormap_algoA_recoToSim
edm::EDGetTokenT< reco::VertexCollection > vertexLabel_algoB
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< reco::BeamSpot > beamSpotLabel
edm::EDGetTokenT< edm::View< reco::Track > > trackLabel_algoA
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
void SetRecoVertex_AlgoB(reco::VertexRef vertex)
Definition: TPtoRecoTrack.h:37
void SetRecoVertex_AlgoA(reco::VertexRef vertex)
Definition: TPtoRecoTrack.h:36
TrackAlgoCompareUtil(const edm::ParameterSet &)
edm::EDGetTokenT< reco::SimToRecoCollection > associatormap_algoA_simToReco
GlobalVector momentum() const
edm::EDGetTokenT< TrackingParticleCollection > trackingParticleLabel_fakes
GlobalPoint position() const
T const * product() const
Definition: Handle.h:81
edm::EDGetTokenT< edm::View< reco::Track > > trackLabel_algoB
edm::EDGetTokenT< reco::TrackToTrackingParticleAssociator > assocLabel_algoB
void SetShared_AlgoB(const float &mB)
Definition: TPtoRecoTrack.h:34
void SetRecoTrack(reco::TrackBaseRef track)
Definition: RecoTracktoTP.h:29
const T & get() const
Definition: EventSetup.h:56
T const * product() const
Definition: ESHandle.h:86
void SetBeamSpot(const math::XYZPoint &bs)
Definition: TPtoRecoTrack.h:39
void SetBeamSpot(const math::XYZPoint &bs)
Definition: RecoTracktoTP.h:31
const Point & position() const
position
Definition: BeamSpot.h:62
void SetRecoTrack_AlgoB(reco::TrackBaseRef track)
Definition: TPtoRecoTrack.h:31
edm::EDGetTokenT< reco::RecoToSimCollection > associatormap_algoB_recoToSim
edm::Ref< TrackingParticleCollection > TrackingParticleRef
Global3DVector GlobalVector
Definition: GlobalVector.h:10
virtual void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override