CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ConversionProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: ConversionProducer
4 // Class: ConversionProducer
5 //
13 //
14 // Original Authors: Hongliang Liu
15 // Created: Thu Mar 13 17:40:48 CDT 2008
16 // $Id: ConversionProducer.cc,v 1.14 2011/12/23 08:13:35 innocent Exp $
17 //
18 //
19 
20 
21 // system include files
22 #include <memory>
23 #include <map>
24 
25 
31 
36 
39 
40 
45 
49 
52 
59 
60 
61 //Kinematic constraint vertex fitter
71 
72 
73 
75  theVertexFinder_(0)
76 
77 {
78  algoName_ = iConfig.getParameter<std::string>( "AlgorithmName" );
79 
80  src_ = iConfig.getParameter<edm::InputTag>("src");
81 
82  maxNumOfTrackInPU_ = iConfig.getParameter<int>("maxNumOfTrackInPU");
83  maxTrackRho_ = iConfig.getParameter<double>("maxTrackRho");
84  maxTrackZ_ = iConfig.getParameter<double>("maxTrackZ");
85 
86  allowTrackBC_ = iConfig.getParameter<bool>("AllowTrackBC");
87  allowD0_ = iConfig.getParameter<bool>("AllowD0");
88  allowDeltaPhi_ = iConfig.getParameter<bool>("AllowDeltaPhi");
89  allowDeltaCot_ = iConfig.getParameter<bool>("AllowDeltaCot");
90  allowMinApproach_ = iConfig.getParameter<bool>("AllowMinApproach");
91  allowOppCharge_ = iConfig.getParameter<bool>("AllowOppCharge");
92 
93  allowVertex_ = iConfig.getParameter<bool>("AllowVertex");
94 
95  bypassPreselGsf_ = iConfig.getParameter<bool>("bypassPreselGsf");
96  bypassPreselEcal_ = iConfig.getParameter<bool>("bypassPreselEcal");
97  bypassPreselEcalEcal_ = iConfig.getParameter<bool>("bypassPreselEcalEcal");
98 
99  deltaEta_ = iConfig.getParameter<double>("deltaEta");
100 
101  halfWayEta_ = iConfig.getParameter<double>("HalfwayEta");//open angle to search track matches with BC
102 
103  d0Cut_ = iConfig.getParameter<double>("d0");
104 
105  usePvtx_ = iConfig.getParameter<bool>("UsePvtx");//if use primary vertices
106 
107  vertexProducer_ = iConfig.getParameter<std::string>("primaryVertexProducer");
108 
109 
110  //Track-cluster matching eta and phi cuts
111  dEtaTkBC_ = iConfig.getParameter<double>("dEtaTrackBC");//TODO research on cut endcap/barrel
112  dPhiTkBC_ = iConfig.getParameter<double>("dPhiTrackBC");
113 
114  bcBarrelCollection_ = iConfig.getParameter<edm::InputTag>("bcBarrelCollection");
115  bcEndcapCollection_ = iConfig.getParameter<edm::InputTag>("bcEndcapCollection");
116 
117  scBarrelProducer_ = iConfig.getParameter<edm::InputTag>("scBarrelProducer");
118  scEndcapProducer_ = iConfig.getParameter<edm::InputTag>("scEndcapProducer");
119 
120  energyBC_ = iConfig.getParameter<double>("EnergyBC");//BC energy threshold
121  energyTotalBC_ = iConfig.getParameter<double>("EnergyTotalBC");//BC energy threshold
122  minSCEt_ = iConfig.getParameter<double>("minSCEt");//super cluster energy threshold
123  dEtacutForSCmatching_ = iConfig.getParameter<double>("dEtacutForSCmatching");// dEta between conversion momentum direction and SC position
124  dPhicutForSCmatching_ = iConfig.getParameter<double>("dPhicutForSCmatching");// dPhi between conversion momentum direction and SC position
125 
126 
127 
128 
129  //Track cuts on left right track: at least one leg reaches ECAL
130  //Left track: must exist, must reach Ecal and match BC, so loose cut on Chi2 and tight on hits
131  //Right track: not necessary to exist (if allowSingleLeg_), not necessary to reach ECAL or match BC, so tight cut on Chi2 and loose on hits
132  maxChi2Left_ = iConfig.getParameter<double>("MaxChi2Left");
133  maxChi2Right_ = iConfig.getParameter<double>("MaxChi2Right");
134  minHitsLeft_ = iConfig.getParameter<int>("MinHitsLeft");
135  minHitsRight_ = iConfig.getParameter<int>("MinHitsRight");
136 
137  //Track Open angle cut on delta cot(theta) and delta phi
138  deltaCotTheta_ = iConfig.getParameter<double>("DeltaCotTheta");
139  deltaPhi_ = iConfig.getParameter<double>("DeltaPhi");
140  minApproachLow_ = iConfig.getParameter<double>("MinApproachLow");
141  minApproachHigh_ = iConfig.getParameter<double>("MinApproachHigh");
142 
143 
144  // if allow single track collection, by default False
145  allowSingleLeg_ = iConfig.getParameter<bool>("AllowSingleLeg");
146  rightBC_ = iConfig.getParameter<bool>("AllowRightBC");
147 
148  //track inner position dz cut, need RECO
149  dzCut_ = iConfig.getParameter<double>("dz");
150  //track analytical cross cut
151  r_cut = iConfig.getParameter<double>("rCut");
152  vtxChi2_ = iConfig.getParameter<double>("vtxChi2");
153 
154 
155  theVertexFinder_ = new ConversionVertexFinder ( iConfig );
156 
157  thettbuilder_ = 0;
158 
159  //output
160  ConvertedPhotonCollection_ = iConfig.getParameter<std::string>("convertedPhotonCollection");
161 
162  produces< reco::ConversionCollection >(ConvertedPhotonCollection_);
163 
164 }
165 
166 
168 {
169 
170  // do anything here that needs to be done at desctruction time
171  // (e.g. close files, deallocate resources etc.)
172  delete theVertexFinder_;
173 }
174 
175 
176 // ------------ method called to produce the data ------------
177 void
179 {
180  using namespace edm;
181 
182  reco::ConversionCollection outputConvPhotonCollection;
183  std::auto_ptr<reco::ConversionCollection> outputConvPhotonCollection_p(new reco::ConversionCollection);
184 
185  //std::cout << " ConversionProducer::produce " << std::endl;
186  //Read multiple track input collections
187 
188  edm::Handle<edm::View<reco::ConversionTrack> > trackCollectionHandle;
189  iEvent.getByLabel(src_,trackCollectionHandle);
190 
191  //build map of ConversionTracks ordered in eta
192  std::multimap<float, edm::Ptr<reco::ConversionTrack> > convTrackMap;
193  for (edm::PtrVector<reco::ConversionTrack>::const_iterator tk_ref = trackCollectionHandle->ptrVector().begin(); tk_ref != trackCollectionHandle->ptrVector().end(); ++tk_ref ){
194  convTrackMap.insert(std::make_pair((*tk_ref)->track()->eta(),*tk_ref));
195  }
196 
199  if (usePvtx_){
200  iEvent.getByLabel(vertexProducer_, vertexHandle);
201  if (!vertexHandle.isValid()) {
202  edm::LogError("ConversionProducer") << "Error! Can't get the product primary Vertex Collection "<< "\n";
203  usePvtx_ = false;
204  }
205  if (usePvtx_)
206  vertexCollection = *(vertexHandle.product());
207  }
208 
209  edm::ESHandle<TransientTrackBuilder> hTransientTrackBuilder;
210  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",hTransientTrackBuilder);
211  thettbuilder_ = hTransientTrackBuilder.product();
212 
213  reco::Vertex the_pvtx;
214  //because the priamry vertex is sorted by quality, the first one is the best
215  if (!vertexCollection.empty())
216  the_pvtx = *(vertexCollection.begin());
217 
218  if (trackCollectionHandle->size()> maxNumOfTrackInPU_){
219  iEvent.put( outputConvPhotonCollection_p, ConvertedPhotonCollection_);
220  return;
221  }
222 
223 
224  // build Super and Basic cluster geometry map to search in eta bounds for clusters
225  std::multimap<double, reco::CaloClusterPtr> basicClusterPtrs;
226  std::multimap<double, reco::CaloClusterPtr> superClusterPtrs;
227 
228 
229  buildSuperAndBasicClusterGeoMap(iEvent,basicClusterPtrs,superClusterPtrs);
230 
231  buildCollection( iEvent, iSetup, convTrackMap, superClusterPtrs, basicClusterPtrs, the_pvtx, outputConvPhotonCollection);//allow empty basicClusterPtrs
232 
233  outputConvPhotonCollection_p->assign(outputConvPhotonCollection.begin(), outputConvPhotonCollection.end());
234  iEvent.put( outputConvPhotonCollection_p, ConvertedPhotonCollection_);
235 
236 }
237 
238 
240  std::multimap<double, reco::CaloClusterPtr>& basicClusterPtrs,
241  std::multimap<double, reco::CaloClusterPtr>& superClusterPtrs){
242 
243  // Get the Super Cluster collection in the Barrel
245  iEvent.getByLabel(scBarrelProducer_,scBarrelHandle);
246  if (!scBarrelHandle.isValid()) {
247  edm::LogError("ConvertedPhotonProducer") << "Error! Can't get the product "<<scBarrelProducer_;
248  }
249 
250  // Get the Super Cluster collection in the Endcap
252  iEvent.getByLabel(scEndcapProducer_,scEndcapHandle);
253  if (!scEndcapHandle.isValid()) {
254  edm::LogError("ConvertedPhotonProducer") << "Error! Can't get the product "<<scEndcapProducer_;
255  }
256 
257 
259  edm::Handle<edm::View<reco::CaloCluster> > bcEndcapHandle;//TODO check cluster type if BasicCluster or PFCluster
260 
261  iEvent.getByLabel( bcBarrelCollection_, bcBarrelHandle);
262  if (!bcBarrelHandle.isValid()) {
263  edm::LogError("ConvertedPhotonProducer") << "Error! Can't get the product "<<bcBarrelCollection_;
264  }
265 
266  iEvent.getByLabel( bcEndcapCollection_, bcEndcapHandle);
267  if (! bcEndcapHandle.isValid()) {
268  edm::LogError("ConvertedPhotonProducer") << "Error! Can't get the product "<<bcEndcapCollection_;
269  }
270 
271  edm::Handle<edm::View<reco::CaloCluster> > bcHandle = bcBarrelHandle;
272  edm::Handle<edm::View<reco::CaloCluster> > scHandle = scBarrelHandle;
273 
274  if ( bcHandle.isValid() ) {
275  for (unsigned jj = 0; jj < 2; ++jj ){
276  for (unsigned ii = 0; ii < bcHandle->size(); ++ii ) {
277  if (bcHandle->ptrAt(ii)->energy()>energyBC_)
278  basicClusterPtrs.insert(std::make_pair(bcHandle->ptrAt(ii)->position().eta(), bcHandle->ptrAt(ii)));
279  }
280  bcHandle = bcEndcapHandle;
281  }
282  }
283 
284 
285  if ( scHandle.isValid() ) {
286  for (unsigned jj = 0; jj < 2; ++jj ){
287  for (unsigned ii = 0; ii < scHandle->size(); ++ii ) {
288  if (scHandle->ptrAt(ii)->energy()>minSCEt_)
289  superClusterPtrs.insert(std::make_pair(scHandle->ptrAt(ii)->position().eta(), scHandle->ptrAt(ii)));
290  }
291  scHandle = scEndcapHandle;
292  }
293  }
294 
295 
296 }
297 
298 
300  const std::multimap<float, edm::Ptr<reco::ConversionTrack> >& allTracks,
301  const std::multimap<double, reco::CaloClusterPtr>& superClusterPtrs,
302  const std::multimap<double, reco::CaloClusterPtr>& basicClusterPtrs,
303  const reco::Vertex& the_pvtx,
304  reco::ConversionCollection & outputConvPhotonCollection){
305 
306 
307  edm::ESHandle<TrackerGeometry> trackerGeomHandle;
308  edm::ESHandle<MagneticField> magFieldHandle;
309  iSetup.get<TrackerDigiGeometryRecord>().get( trackerGeomHandle );
310  iSetup.get<IdealMagneticFieldRecord>().get( magFieldHandle );
311 
312  const TrackerGeometry* trackerGeom = trackerGeomHandle.product();
313  const MagneticField* magField = magFieldHandle.product();
314 
315 // std::vector<math::XYZPointF> trackImpactPosition;
316 // trackImpactPosition.reserve(allTracks.size());//track impact position at ECAL
317 // std::vector<bool> trackValidECAL;//Does this track reach ECAL basic cluster (reach ECAL && match with BC)
318 // trackValidECAL.assign(allTracks.size(), false);
319 //
320 // std::vector<reco::CaloClusterPtr> trackMatchedBC;
321 // reco::CaloClusterPtr empty_bc;
322 // trackMatchedBC.assign(allTracks.size(), empty_bc);//TODO find a better way to avoid copy constructor
323 //
324 // std::vector<int> bcHandleId;//the associated BC handle id, -1 invalid, 0 barrel 1 endcap
325 // bcHandleId.assign(allTracks.size(), -1);
326 
327  // not used std::multimap<double, int> trackInnerEta;//Track innermost state Eta map to TrackRef index, to be used in track pair sorting
328 
329  std::map<edm::Ptr<reco::ConversionTrack>, math::XYZPointF> trackImpactPosition;
330  std::map<edm::Ptr<reco::ConversionTrack>, reco::CaloClusterPtr> trackMatchedBC;
331 
332  ConversionHitChecker hitChecker;
333 
334 
335  //2 propagate all tracks into ECAL, record its eta and phi
336 
337  for (std::multimap<float, edm::Ptr<reco::ConversionTrack> >::const_iterator tk_ref = allTracks.begin(); tk_ref != allTracks.end(); ++tk_ref ){
338  const reco::Track* tk = tk_ref->second->trackRef().get() ;
339 
340 
341  //check impact position then match with BC
342  math::XYZPointF ew;
343  if ( getTrackImpactPosition(tk, trackerGeom, magField, ew) ){
344  trackImpactPosition[tk_ref->second] = ew;
345 
346  reco::CaloClusterPtr closest_bc;//the closest matching BC to track
347 
348  if ( getMatchedBC(basicClusterPtrs, ew, closest_bc) ){
349  trackMatchedBC[tk_ref->second] = closest_bc;
350  }
351  }
352  }
353 
354 
355 
356  //3. pair up tracks:
357  //TODO it is k-Closest pair of point problem
358  //std::cout << " allTracks.size() " << allTracks.size() << std::endl;
359  for(std::multimap<float, edm::Ptr<reco::ConversionTrack> >::const_iterator ll = allTracks.begin(); ll != allTracks.end(); ++ll ) {
360  bool track1HighPurity=true;
361  //std::cout << " Loop on allTracks " << std::endl;
362  const edm::RefToBase<reco::Track> & left = ll->second->trackRef();
363 
364 
365  //TODO: This is a workaround, should be fixed with a proper function in the TTBuilder
366  //(Note that the TrackRef and GsfTrackRef versions of the constructor are needed
367  // to properly get refit tracks in the output vertex)
368  reco::TransientTrack ttk_l;
369  if (dynamic_cast<const reco::GsfTrack*>(left.get())) {
370  ttk_l = thettbuilder_->build(left.castTo<reco::GsfTrackRef>());
371  }
372  else {
373  ttk_l = thettbuilder_->build(left.castTo<reco::TrackRef>());
374  }
375 
377  // if ((allowTrackBC_ && !trackValidECAL[ll-allTracks.begin()]) )//this Left leg should have valid BC
378  // continue;
379 
380 
381  if (the_pvtx.isValid()){
382  if (!(trackD0Cut(left, the_pvtx))) track1HighPurity=false;
383  } else {
384  if (!(trackD0Cut(left))) track1HighPurity=false;
385  }
386 
387  std::vector<int> right_candidates;//store all right legs passed the cut (theta/approach and ref pair)
388  std::vector<double> right_candidate_theta, right_candidate_approach;
389  std::vector<std::pair<bool, reco::Vertex> > vertex_candidates;
390 
391  //inner loop only over tracks between eta and eta + deltaEta of the first track
392  float etasearch = ll->first + deltaEta_;
393  std::multimap<float, edm::Ptr<reco::ConversionTrack> >::const_iterator rr = ll;
394  ++rr;
395  for (; rr != allTracks.lower_bound(etasearch); ++rr ) {
396  bool track2HighPurity = true;
397  bool highPurityPair = true;
398 
399  const edm::RefToBase<reco::Track> & right = rr->second->trackRef();
400 
401 
402  //TODO: This is a workaround, should be fixed with a proper function in the TTBuilder
403  reco::TransientTrack ttk_r;
404  if (dynamic_cast<const reco::GsfTrack*>(right.get())) {
405  ttk_r = thettbuilder_->build(right.castTo<reco::GsfTrackRef>());
406  }
407  else {
408  ttk_r = thettbuilder_->build(right.castTo<reco::TrackRef>());
409  }
410  //std::cout << " This track is " << right->algoName() << std::endl;
411 
412 
413  //all vertexing preselection should go here
414 
415  //check for opposite charge
416  if ( allowOppCharge_ && (left->charge()*right->charge() > 0) )
417  continue; //same sign, reject pair
418 
420  //if ( (allowTrackBC_ && !trackValidECAL[rr-allTracks.begin()] && rightBC_) )// if right track matches ECAL
421  // continue;
422 
423 
424  double approachDist = -999.;
425  //apply preselection to track pair, overriding preselection for gsf+X or ecalseeded+X pairs if so configured
426  bool preselected = preselectTrackPair(ttk_l,ttk_r, approachDist);
427  preselected = preselected || (bypassPreselGsf_ && (left->algo()==29 || right->algo()==29));
428  preselected = preselected || (bypassPreselEcal_ && (left->algo()==15 || right->algo()==15 || left->algo()==16 || right->algo()==16));
429  preselected = preselected || (bypassPreselEcalEcal_ && (left->algo()==15 || left->algo()==16) && (right->algo()==15 || right->algo()==16));
430 
431  if (!preselected) {
432  continue;
433  }
434 
435  //do the actual vertex fit
436  reco::Vertex theConversionVertex;//by default it is invalid
437  bool goodVertex = checkVertex(ttk_l, ttk_r, magField, theConversionVertex);
438 
439  //bail as early as possible in case the fit didn't return a good vertex
440  if (!goodVertex) {
441  continue;
442  }
443 
444 
445 
446  //track pair pass the quality cut
447  if ( !( (trackQualityFilter(left, true) && trackQualityFilter(right, false))
448  || (trackQualityFilter(left, false) && trackQualityFilter(right, true)) ) ) {
449  highPurityPair=false;
450  }
451 
452  if (the_pvtx.isValid()){
453  if (!(trackD0Cut(right, the_pvtx))) track2HighPurity=false;
454  } else {
455  if (!(trackD0Cut(right))) track2HighPurity=false;
456  }
457 
458 
459  //if all cuts passed, go ahead to make conversion candidates
460  std::vector<edm::RefToBase<reco::Track> > trackPairRef;
461  trackPairRef.push_back(left);//left track
462  trackPairRef.push_back(right);//right track
463 
464  std::vector<math::XYZVectorF> trackPin;
465  std::vector<math::XYZVectorF> trackPout;
466  std::vector<math::XYZPointF> trackInnPos;
467  std::vector<uint8_t> nHitsBeforeVtx;
468  std::vector<Measurement1DFloat> dlClosestHitToVtx;
469 
470  if (left->extra().isNonnull() && right->extra().isNonnull()){//only available on TrackExtra
471  trackInnPos.push_back( toFConverterP(left->innerPosition()));
472  trackInnPos.push_back( toFConverterP(right->innerPosition()));
473  trackPin.push_back(toFConverterV(left->innerMomentum()));
474  trackPin.push_back(toFConverterV(right->innerMomentum()));
475  trackPout.push_back(toFConverterV(left->outerMomentum()));
476  trackPout.push_back(toFConverterV(right->outerMomentum()));
477  }
478 
479  if (ll->second->trajRef().isNonnull() && rr->second->trajRef().isNonnull()) {
480  std::pair<uint8_t,Measurement1DFloat> leftWrongHits = hitChecker.nHitsBeforeVtx(*ll->second->trajRef().get(),theConversionVertex);
481  std::pair<uint8_t,Measurement1DFloat> rightWrongHits = hitChecker.nHitsBeforeVtx(*rr->second->trajRef().get(),theConversionVertex);
482  nHitsBeforeVtx.push_back(leftWrongHits.first);
483  nHitsBeforeVtx.push_back(rightWrongHits.first);
484  dlClosestHitToVtx.push_back(leftWrongHits.second);
485  dlClosestHitToVtx.push_back(rightWrongHits.second);
486  }
487 
488  uint8_t nSharedHits = hitChecker.nSharedHits(*left.get(),*right.get());
489 
490 
491  //if using kinematic fit, check with chi2 post cut
492  if (theConversionVertex.isValid()){
493  const float chi2Prob = ChiSquaredProbability(theConversionVertex.chi2(), theConversionVertex.ndof());
494  if (chi2Prob<vtxChi2_) highPurityPair=false;
495  }
496 
497  //std::cout << " highPurityPair after vertex cut " << highPurityPair << std::endl;
498  std::vector<math::XYZPointF> trkPositionAtEcal;
499  std::vector<reco::CaloClusterPtr> matchingBC;
500 
501  if (allowTrackBC_){//TODO find out the BC ptrs if not doing matching, otherwise, leave it empty
502  //const int lbc_handle = bcHandleId[ll-allTracks.begin()],
503  // rbc_handle = bcHandleId[rr-allTracks.begin()];
504 
505  std::map<edm::Ptr<reco::ConversionTrack>, math::XYZPointF>::const_iterator trackImpactPositionLeft = trackImpactPosition.find(ll->second);
506  std::map<edm::Ptr<reco::ConversionTrack>, math::XYZPointF>::const_iterator trackImpactPositionRight = trackImpactPosition.find(rr->second);
507  std::map<edm::Ptr<reco::ConversionTrack>, reco::CaloClusterPtr>::const_iterator trackMatchedBCLeft = trackMatchedBC.find(ll->second);
508  std::map<edm::Ptr<reco::ConversionTrack>, reco::CaloClusterPtr>::const_iterator trackMatchedBCRight = trackMatchedBC.find(rr->second);
509 
510  if (trackImpactPositionLeft!=trackImpactPosition.end()) {
511  trkPositionAtEcal.push_back(trackImpactPositionLeft->second);//left track
512  }
513  else {
514  trkPositionAtEcal.push_back(math::XYZPointF());//left track
515  }
516  if (trackImpactPositionRight!=trackImpactPosition.end()) {//second track ECAL position may be invalid
517  trkPositionAtEcal.push_back(trackImpactPositionRight->second);
518  }
519 
520  double total_e_bc = 0.;
521  if (trackMatchedBCLeft!=trackMatchedBC.end()) {
522  matchingBC.push_back(trackMatchedBCLeft->second);//left track
523  total_e_bc += trackMatchedBCLeft->second->energy();
524  }
525  else {
526  matchingBC.push_back( reco::CaloClusterPtr() );//left track
527  }
528  if (trackMatchedBCRight!=trackMatchedBC.end()) {//second track ECAL position may be invalid
529  matchingBC.push_back(trackMatchedBCRight->second);
530  total_e_bc += trackMatchedBCRight->second->energy();
531  }
532 
533  if (total_e_bc<energyTotalBC_) {
534  highPurityPair = false;
535  }
536 
537 
538  }
539  //signature cuts, then check if vertex, then post-selection cuts
540  highPurityPair = highPurityPair && track1HighPurity && track2HighPurity && goodVertex && checkPhi(left, right, trackerGeom, magField, theConversionVertex) ;
541 
542 
544  /*
545  for ( std::vector<edm::RefToBase<reco::Track> >::iterator iTk=trackPairRef.begin(); iTk!=trackPairRef.end(); iTk++) {
546  math::XYZPointF impPos;
547  if ( getTrackImpactPosition(*iTk, trackerGeom, magField, impPos) ) {
548 
549  }
550 
551  }
552  */
553 
554  const float minAppDist = approachDist;
556  float dummy=0;
558  reco::Conversion newCandidate(scPtrVec, trackPairRef, trkPositionAtEcal, theConversionVertex, matchingBC, minAppDist, trackInnPos, trackPin, trackPout, nHitsBeforeVtx, dlClosestHitToVtx, nSharedHits, dummy, algo );
559  // Fill in scPtrVec with the macthing SC
560  if ( matchingSC ( superClusterPtrs, newCandidate, scPtrVec) )
561  newCandidate.setMatchingSuperCluster( scPtrVec);
562 
563  //std::cout << " ConversionProducer scPtrVec.size " << scPtrVec.size() << std::endl;
564 
565  newCandidate.setQuality(reco::Conversion::highPurity, highPurityPair);
566  bool generalTracksOnly = ll->second->isTrackerOnly() && rr->second->isTrackerOnly() && !dynamic_cast<const reco::GsfTrack*>(ll->second->trackRef().get()) && !dynamic_cast<const reco::GsfTrack*>(rr->second->trackRef().get());
567  bool arbitratedEcalSeeded = ll->second->isArbitratedEcalSeeded() && rr->second->isArbitratedEcalSeeded();
568  bool arbitratedMerged = ll->second->isArbitratedMerged() && rr->second->isArbitratedMerged();
569  bool arbitratedMergedEcalGeneral = ll->second->isArbitratedMergedEcalGeneral() && rr->second->isArbitratedMergedEcalGeneral();
570 
571  newCandidate.setQuality(reco::Conversion::generalTracksOnly, generalTracksOnly);
572  newCandidate.setQuality(reco::Conversion::arbitratedEcalSeeded, arbitratedEcalSeeded);
573  newCandidate.setQuality(reco::Conversion::arbitratedMerged, arbitratedMerged);
574  newCandidate.setQuality(reco::Conversion::arbitratedMergedEcalGeneral, arbitratedMergedEcalGeneral);
575 
576  outputConvPhotonCollection.push_back(newCandidate);
577 
578  }
579 
580  }
581 
582 
583 
584 
585 
586 
587 }
588 
589 
590 
591 
592 
593 //
594 // member functions
595 //
596 
598  bool pass = true;
599  if (isLeft){
600  pass = (ref->normalizedChi2() < maxChi2Left_ && ref->found() >= minHitsLeft_);
601  } else {
602  pass = (ref->normalizedChi2() < maxChi2Right_ && ref->found() >= minHitsRight_);
603  }
604 
605  return pass;
606 }
607 
609  //NOTE if not allow d0 cut, always true
610  return ((!allowD0_) || !(ref->d0()*ref->charge()/ref->d0Error()<d0Cut_));
611 }
612 
614  //
615  return ((!allowD0_) || !(-ref->dxy(the_pvtx.position())*ref->charge()/ref->dxyError()<d0Cut_));
616 }
617 
618 
620  const TrackerGeometry* trackerGeom, const MagneticField* magField,
621  math::XYZPointF& ew){
622 
623  PropagatorWithMaterial propag( alongMomentum, 0.000511, magField );
624 
626  new BoundCylinder( GlobalPoint(0.,0.,0.), TkRotation<float>(),
627  SimpleCylinderBounds( 129, 129, -320.5, 320.5 ) ) );
628  const float epsilon = 0.001;
629  Surface::RotationType rot; // unit rotation matrix
630  const float barrelRadius = 129.f;
631  const float barrelHalfLength = 270.9f;
632  const float endcapRadius = 171.1f;
633  const float endcapZ = 320.5f;
635  SimpleCylinderBounds( barrelRadius-epsilon, barrelRadius+epsilon, -barrelHalfLength, barrelHalfLength)));
636  ReferenceCountingPointer<BoundDisk> theNegativeEtaEndcap_(
637  new BoundDisk( Surface::PositionType( 0, 0, -endcapZ), rot,
638  SimpleDiskBounds( 0, endcapRadius, -epsilon, epsilon)));
639  ReferenceCountingPointer<BoundDisk> thePositiveEtaEndcap_(
640  new BoundDisk( Surface::PositionType( 0, 0, endcapZ), rot,
641  SimpleDiskBounds( 0, endcapRadius, -epsilon, epsilon)));
642 
643  //const TrajectoryStateOnSurface myTSOS = trajectoryStateTransform::innerStateOnSurface(*(*ref), *trackerGeom, magField);
644  const TrajectoryStateOnSurface myTSOS = trajectoryStateTransform::outerStateOnSurface(*tk_ref, *trackerGeom, magField);
645  TrajectoryStateOnSurface stateAtECAL;
646  stateAtECAL = propag.propagate(myTSOS, *theBarrel_);
647  if (!stateAtECAL.isValid() || ( stateAtECAL.isValid() && fabs(stateAtECAL.globalPosition().eta() ) >1.479 ) ) {
648  //endcap propagator
649  if (myTSOS.globalPosition().eta() > 0.) {
650  stateAtECAL = propag.propagate(myTSOS, *thePositiveEtaEndcap_);
651  } else {
652  stateAtECAL = propag.propagate(myTSOS, *theNegativeEtaEndcap_);
653  }
654  }
655  if (stateAtECAL.isValid()){
656  ew = stateAtECAL.globalPosition();
657  return true;
658  }
659  else
660  return false;
661 }
662 
663 
664 
665 
666 bool ConversionProducer::matchingSC(const std::multimap<double, reco::CaloClusterPtr>& scMap,
667  reco::Conversion& aConv,
668  // reco::CaloClusterPtr& mSC){
670 
671  // double dRMin=999.;
672  double detaMin=999.;
673  double dphiMin=999.;
675  for (std::multimap<double, reco::CaloClusterPtr>::const_iterator scItr = scMap.begin(); scItr != scMap.end(); scItr++) {
676  const reco::CaloClusterPtr& sc = scItr->second;
677  const double delta_phi = reco::deltaPhi( aConv.refittedPairMomentum().phi(), sc->phi());
678  double sceta = sc->eta();
679  double conveta = etaTransformation(aConv.refittedPairMomentum().eta(), aConv.zOfPrimaryVertexFromTracks() );
680  const double delta_eta = fabs(conveta - sceta);
681  if ( fabs(delta_eta) < fabs(detaMin) && fabs(delta_phi) < fabs(dphiMin) ) {
682  detaMin= fabs(delta_eta);
683  dphiMin= fabs(delta_phi);
684  match=sc;
685  }
686  }
687 
688  if ( fabs(detaMin) < dEtacutForSCmatching_ && fabs(dphiMin) < dPhicutForSCmatching_ ) {
689  mSC.push_back(match);
690  return true;
691  } else
692  return false;
693 }
694 
695 bool ConversionProducer::getMatchedBC(const std::multimap<double, reco::CaloClusterPtr>& bcMap,
696  const math::XYZPointF& trackImpactPosition,
697  reco::CaloClusterPtr& closestBC){
698  const double track_eta = trackImpactPosition.eta();
699  const double track_phi = trackImpactPosition.phi();
700 
701  double min_eta = 999., min_phi = 999.;
702  reco::CaloClusterPtr closest_bc;
703  for (std::multimap<double, reco::CaloClusterPtr>::const_iterator bc = bcMap.lower_bound(track_eta - halfWayEta_);
704  bc != bcMap.upper_bound(track_eta + halfWayEta_); ++bc){//use eta map to select possible BC collection then loop in
705  const reco::CaloClusterPtr& ebc = bc->second;
706  const double delta_eta = track_eta-(ebc->position().eta());
707  const double delta_phi = reco::deltaPhi(track_phi, (ebc->position().phi()));
708  if (fabs(delta_eta)<dEtaTkBC_ && fabs(delta_phi)<dPhiTkBC_){
709  if (fabs(min_eta)>fabs(delta_eta) && fabs(min_phi)>fabs(delta_phi)){//take the closest to track BC
710  min_eta = delta_eta;
711  min_phi = delta_phi;
712  closest_bc = bc->second;
713  //TODO check if min_eta>delta_eta but min_phi<delta_phi
714  }
715  }
716  }
717 
718  if (min_eta < 999.){
719  closestBC = closest_bc;
720  return true;
721  } else
722  return false;
723 }
724 
725 
726 
727 
728 
729 
730 //check track open angle of phi at vertex
732  const TrackerGeometry* trackerGeom, const MagneticField* magField,
733  const reco::Vertex& vtx){
734  if (!allowDeltaPhi_)
735  return true;
736  //if track has innermost momentum, check with innermost phi
737  //if track also has valid vertex, propagate to vertex then calculate phi there
738  //if track has no innermost momentum, just return true, because track->phi() makes no sense
739  if (tk_l->extra().isNonnull() && tk_r->extra().isNonnull()){
740  double iphi1 = tk_l->innerMomentum().phi(), iphi2 = tk_r->innerMomentum().phi();
741  if (vtx.isValid()){
742  PropagatorWithMaterial propag( anyDirection, 0.000511, magField );
743 
744  double recoPhoR = vtx.position().Rho();
747  SimpleCylinderBounds( recoPhoR-0.001, recoPhoR+0.001, -fabs(vtx.position().z()), fabs(vtx.position().z()))));
749  new BoundDisk( Surface::PositionType( 0, 0, vtx.position().z()), rot,
750  SimpleDiskBounds( 0, recoPhoR, -0.001, 0.001)));
751 
752  const TrajectoryStateOnSurface myTSOS1 = trajectoryStateTransform::innerStateOnSurface(*tk_l, *trackerGeom, magField);
753  const TrajectoryStateOnSurface myTSOS2 = trajectoryStateTransform::innerStateOnSurface(*tk_r, *trackerGeom, magField);
754  TrajectoryStateOnSurface stateAtVtx1, stateAtVtx2;
755  stateAtVtx1 = propag.propagate(myTSOS1, *theBarrel_);
756  if (!stateAtVtx1.isValid() ) {
757  stateAtVtx1 = propag.propagate(myTSOS1, *theDisk_);
758  }
759  if (stateAtVtx1.isValid()){
760  iphi1 = stateAtVtx1.globalDirection().phi();
761  }
762  stateAtVtx2 = propag.propagate(myTSOS2, *theBarrel_);
763  if (!stateAtVtx2.isValid() ) {
764  stateAtVtx2 = propag.propagate(myTSOS2, *theDisk_);
765  }
766  if (stateAtVtx2.isValid()){
767  iphi2 = stateAtVtx2.globalDirection().phi();
768  }
769  }
770  const double dPhi = reco::deltaPhi(iphi1, iphi2);
771  return (fabs(dPhi) < deltaPhi_);
772  } else {
773  return true;
774  }
775 }
776 
778  double& appDist) {
779 
780 
781  double dCotTheta = 1./tan(ttk_l.track().innerMomentum().theta()) - 1./tan(ttk_r.track().innerMomentum().theta());
782  if (allowDeltaCot_ && (std::abs(dCotTheta) > deltaCotTheta_)) {
783  return false;
784  }
785 
786  //non-conversion hypothesis, reject prompt track pairs
787  ClosestApproachInRPhi closest;
789  if (!closest.status()) {
790  return false;
791  }
792 
793  if (closest.crossingPoint().perp() < r_cut) {
794  return false;
795  }
796 
797 
798  //compute tangent point btw tracks (conversion hypothesis)
799  TangentApproachInRPhi tangent;
801  if (!tangent.status()) {
802  return false;
803  }
804 
805  GlobalPoint tangentPoint = tangent.crossingPoint();
806  double rho = tangentPoint.perp();
807 
808  //reject candidates well outside of tracker bounds
809  if (rho > maxTrackRho_) {
810  return false;
811  }
812 
813  if (std::abs(tangentPoint.z()) > maxTrackZ_) {
814  return false;
815  }
816 
817  std::pair<GlobalTrajectoryParameters,GlobalTrajectoryParameters> trajs = tangent.trajectoryParameters();
818 
819  //very large separation in z, no hope
820  if (std::abs(trajs.first.position().z() - trajs.second.position().z()) > dzCut_) {
821  return false;
822  }
823 
824 
825  float minApproach = tangent.perpdist();
826  appDist = minApproach;
827 
828  if (allowMinApproach_ && (minApproach < minApproachLow_ || minApproach > minApproachHigh_) ) {
829  return false;
830  }
831 
832  return true;
833 
834 
835 }
836 
839 
840  const reco::CaloClusterPtr& bc_l = ll.second;//can be null, so check isNonnull()
841  const reco::CaloClusterPtr& bc_r = rr.second;
842 
843  //The cuts should be ordered by considering if takes time and if cuts off many fakes
844  if (allowTrackBC_){
845  //check energy of BC
846  double total_e_bc = 0;
847  if (bc_l.isNonnull()) total_e_bc += bc_l->energy();
848  if (rightBC_)
849  if (bc_r.isNonnull())
850  total_e_bc += bc_r->energy();
851 
852  if (total_e_bc < energyTotalBC_) return false;
853  }
854 
855  return true;
856 }
857 
858 
859 
860 //because reco::vertex uses track ref, so have to keep them
862  const MagneticField* magField,
863  reco::Vertex& the_vertex){
864  bool found = false;
865 
866  std::vector<reco::TransientTrack> pair;
867  pair.push_back(ttk_l);
868  pair.push_back(ttk_r);
869 
870  found = theVertexFinder_->run(pair, the_vertex);
871 
872 
873 
874  return found;
875 }
876 
877 
878 
879 double ConversionProducer::etaTransformation( float EtaParticle , float Zvertex) {
880 
881  //---Definitions
882  const float PI = 3.1415927;
883 
884  //---Definitions for ECAL
885  const float R_ECAL = 136.5;
886  const float Z_Endcap = 328.0;
887  const float etaBarrelEndcap = 1.479;
888 
889  //---ETA correction
890 
891  float Theta = 0.0 ;
892  float ZEcal = R_ECAL*sinh(EtaParticle)+Zvertex;
893 
894  if(ZEcal != 0.0) Theta = atan(R_ECAL/ZEcal);
895  if(Theta<0.0) Theta = Theta+PI ;
896  double ETA = - log(tan(0.5*Theta));
897 
898  if( fabs(ETA) > etaBarrelEndcap )
899  {
900  float Zend = Z_Endcap ;
901  if(EtaParticle<0.0 ) Zend = -Zend ;
902  float Zlen = Zend - Zvertex ;
903  float RR = Zlen/sinh(EtaParticle);
904  Theta = atan(RR/Zend);
905  if(Theta<0.0) Theta = Theta+PI ;
906  ETA = - log(tan(0.5*Theta));
907  }
908  //---Return the result
909  return ETA;
910  //---end
911 }
912 
T getParameter(std::string const &) const
unsigned int maxNumOfTrackInPU_
ConversionVertexFinder * theVertexFinder_
double d0Error() const
error on d0
Definition: TrackBase.h:211
double delta_eta(double eta1, double eta2)
Definition: AnglesUtil.h:98
edm::InputTag scBarrelProducer_
T perp() const
Definition: PV3DBase.h:71
double d0() const
dxy parameter in perigee convention (d0 = - dxy)
Definition: TrackBase.h:123
bool getTrackImpactPosition(const reco::Track *tk_ref, const TrackerGeometry *trackerGeom, const MagneticField *magField, math::XYZPointF &ew)
const TrackExtraRef & extra() const
reference to &quot;extra&quot; object
Definition: Track.h:97
TrajectoryStateOnSurface innerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field)
std::pair< GlobalTrajectoryParameters, GlobalTrajectoryParameters > trajectoryParameters() const
#define PI
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:111
edm::InputTag bcBarrelCollection_
double dxyError() const
error on dxy
Definition: TrackBase.h:209
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:137
std::pair< uint8_t, Measurement1DFloat > nHitsBeforeVtx(const Trajectory &traj, const reco::Vertex &vtx, double sigmaTolerance=3.0) const
const TransientTrackBuilder * thettbuilder_
void setQuality(ConversionQuality q, bool b)
Definition: Conversion.h:253
bool isValid() const
Tells whether the vertex is valid.
Definition: Vertex.h:61
Definition: DDAxes.h:10
Geom::Phi< T > phi() const
Definition: PV3DBase.h:68
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
void buildCollection(edm::Event &iEvent, const edm::EventSetup &iSetup, const std::multimap< float, edm::Ptr< reco::ConversionTrack > > &allTracks, const std::multimap< double, reco::CaloClusterPtr > &superClusterPtrs, const std::multimap< double, reco::CaloClusterPtr > &basicClusterPtrs, const reco::Vertex &the_pvtx, reco::ConversionCollection &outputConvPhotonCollection)
reco::TransientTrack build(const reco::Track *p) const
double zOfPrimaryVertexFromTracks(const math::XYZPoint &myBeamSpot=math::XYZPoint()) const
Definition: Conversion.h:149
#define abs(x)
Definition: mlp_lapack.h:159
static ConversionAlgorithm algoByName(const std::string &name)
Definition: Conversion.cc:167
void buildSuperAndBasicClusterGeoMap(const edm::Event &, std::multimap< double, reco::CaloClusterPtr > &basicClusterPtrs, std::multimap< double, reco::CaloClusterPtr > &superClusterPtrs)
GlobalPoint globalPosition() const
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float > > XYZPointF
point in space with cartesian internal representation
Definition: Point3D.h:11
const Point & position() const
position
Definition: Vertex.h:93
tuple vertexCollection
math::XYZVectorF refittedPairMomentum() const
Conversion tracks momentum from the tracks refitted with vertex constraint.
Definition: Conversion.cc:262
std::vector< Conversion > ConversionCollection
collectin of Conversion objects
Definition: ConversionFwd.h:9
edm::InputTag scEndcapProducer_
TrajectoryStateOnSurface innermostMeasurementState() const
bool matchingSC(const std::multimap< double, reco::CaloClusterPtr > &scMap, reco::Conversion &conv, reco::CaloClusterPtrVector &mSC)
virtual GlobalPoint crossingPoint() const
virtual GlobalPoint crossingPoint() const
const math::XYZPoint & innerPosition() const
position of the innermost hit
Definition: Track.h:42
TrackAlgorithm algo() const
Definition: TrackBase.h:332
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:152
#define ETA
int iEvent
Definition: GenABIO.cc:243
double dPhi(double phi1, double phi2)
Definition: JetUtil.h:30
TransientVertex run(std::vector< reco::TransientTrack > pair)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
math::XYZVectorF toFConverterV(const math::XYZVector &val)
T z() const
Definition: PV3DBase.h:63
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
double chi2() const
chi-squares
Definition: Vertex.h:82
float ChiSquaredProbability(double chiSquared, double nrDOF)
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
uint8_t nSharedHits(const reco::Track &trk1, const reco::Track &trk2) const
bool getMatchedBC(const std::multimap< double, reco::CaloClusterPtr > &bcMap, const math::XYZPointF &trackImpactPosition, reco::CaloClusterPtr &closestBC)
double ndof() const
Definition: Vertex.h:89
double delta_phi(double ph11, double phi2)
Definition: AnglesUtil.h:91
bool checkTrackPair(const std::pair< edm::RefToBase< reco::Track >, reco::CaloClusterPtr > &ll, const std::pair< edm::RefToBase< reco::Track >, reco::CaloClusterPtr > &rr)
bool trackQualityFilter(const edm::RefToBase< reco::Track > &ref, bool isLeft)
virtual TrajectoryStateOnSurface propagate(const TrajectoryStateOnSurface &tsos, const Plane &plane) const
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:12
const math::XYZVector & outerMomentum() const
momentum vector at the outermost hit position
Definition: Track.h:49
virtual bool calculate(const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb)
const Track & track() const
bool trackD0Cut(const edm::RefToBase< reco::Track > &ref)
const T & get() const
Definition: EventSetup.h:55
REF castTo() const
cast to a concrete type
Definition: RefToBase.h:241
TrajectoryStateOnSurface outerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field)
T const * product() const
Definition: ESHandle.h:62
std::string ConvertedPhotonCollection_
static const float etaBarrelEndcap
static const float Z_Endcap
T eta() const
Definition: PV3DBase.h:75
virtual bool status() const
edm::InputTag bcEndcapCollection_
const math::XYZVector & innerMomentum() const
momentum vector at the innermost hit position
Definition: Track.h:45
unsigned short found() const
Number of valid hits on track.
Definition: Track.h:100
static const float R_ECAL
double etaTransformation(float EtaParticle, float Zvertex)
int charge() const
track electric charge
Definition: TrackBase.h:113
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:6
math::XYZPointF toFConverterP(const math::XYZPoint &val)
LimitAlgo * algo
Definition: Combine.cc:60
const double epsilon
virtual void produce(edm::Event &, const edm::EventSetup &)
bool preselectTrackPair(const reco::TransientTrack &ttk_l, const reco::TransientTrack &ttk_r, double &appDist)
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:121
virtual bool calculate(const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb)
value_type const * get() const
Definition: RefToBase.h:212
void setMatchingSuperCluster(const reco::CaloClusterPtrVector &sc)
Definition: Conversion.h:177
bool checkPhi(const edm::RefToBase< reco::Track > &tk_l, const edm::RefToBase< reco::Track > &tk_r, const TrackerGeometry *trackerGeom, const MagneticField *magField, const reco::Vertex &the_vertex)
virtual bool status() const
GlobalVector globalDirection() const
ConversionProducer(const edm::ParameterSet &)
bool checkVertex(const reco::TransientTrack &ttk_l, const reco::TransientTrack &ttk_r, const MagneticField *magField, reco::Vertex &the_vertex)