CMS 3D CMS Logo

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

#include <ConversionProducer.h>

Inheritance diagram for ConversionProducer:
edm::stream::EDProducer<>

Public Member Functions

 ConversionProducer (const edm::ParameterSet &)
 
 ~ConversionProducer () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndRuns () const final
 

Private Types

typedef math::XYZPointF Point
 
typedef std::vector< PointPointCollection
 

Private Member Functions

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)
 
void buildSuperAndBasicClusterGeoMap (const edm::Event &, std::multimap< double, reco::CaloClusterPtr > &basicClusterPtrs, std::multimap< double, reco::CaloClusterPtr > &superClusterPtrs)
 
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)
 
bool checkTrackPair (const std::pair< edm::RefToBase< reco::Track >, reco::CaloClusterPtr > &ll, const std::pair< edm::RefToBase< reco::Track >, reco::CaloClusterPtr > &rr)
 
bool checkVertex (const reco::TransientTrack &ttk_l, const reco::TransientTrack &ttk_r, const MagneticField *magField, reco::Vertex &the_vertex)
 
double etaTransformation (float EtaParticle, float Zvertex)
 
bool getMatchedBC (const std::multimap< double, reco::CaloClusterPtr > &bcMap, const math::XYZPointF &trackImpactPosition, reco::CaloClusterPtr &closestBC)
 
bool getTrackImpactPosition (const reco::Track *tk_ref, const TrackerGeometry *trackerGeom, const MagneticField *magField, math::XYZPointF &ew)
 
bool matchingSC (const std::multimap< double, reco::CaloClusterPtr > &scMap, reco::Conversion &conv, reco::CaloClusterPtrVector &mSC)
 
bool preselectTrackPair (const reco::TransientTrack &ttk_l, const reco::TransientTrack &ttk_r, double &appDist)
 
void produce (edm::Event &, const edm::EventSetup &) override
 
math::XYZPointF toFConverterP (const math::XYZPoint &val)
 
math::XYZVectorF toFConverterV (const math::XYZVector &val)
 
bool trackD0Cut (const edm::RefToBase< reco::Track > &ref)
 
bool trackD0Cut (const edm::RefToBase< reco::Track > &ref, const reco::Vertex &the_pvtx)
 
bool trackQualityFilter (const edm::RefToBase< reco::Track > &ref, bool isLeft)
 

Private Attributes

std::string algoName_
 
bool allowD0_
 
bool allowDeltaCot_
 
bool allowDeltaPhi_
 
bool allowMinApproach_
 
bool allowOppCharge_
 
bool allowSingleLeg_
 
bool allowTrackBC_
 
bool allowVertex_
 
edm::EDGetTokenT< edm::View< reco::CaloCluster > > bcBarrelCollection_
 
edm::EDGetTokenT< edm::View< reco::CaloCluster > > bcEndcapCollection_
 
bool bypassPreselEcal_
 
bool bypassPreselEcalEcal_
 
bool bypassPreselGsf_
 
std::string ConvertedPhotonCollection_
 
double d0Cut_
 
double deltaCotTheta_
 
double deltaEta_
 
double deltaPhi_
 
double dEtacutForSCmatching_
 
double dEtaTkBC_
 
double dPhicutForSCmatching_
 
double dPhiTkBC_
 
double dzCut_
 
double energyBC_
 
double energyTotalBC_
 
double halfWayEta_
 
double halfWayPhi_
 
double maxChi2Left_
 
double maxChi2Right_
 
unsigned int maxNumOfTrackInPU_
 
double maxTrackRho_
 
double maxTrackZ_
 
double minApproachHigh_
 
double minApproachLow_
 
double minHitsLeft_
 
double minHitsRight_
 
double minSCEt_
 
double r_cut
 
bool rightBC_
 
edm::EDGetTokenT< edm::View< reco::CaloCluster > > scBarrelProducer_
 
edm::EDGetTokenT< edm::View< reco::CaloCluster > > scEndcapProducer_
 
edm::EDGetTokenT< edm::View< reco::ConversionTrack > > src_
 
const TransientTrackBuilderthettbuilder_
 
ConversionVertexFindertheVertexFinder_
 
bool usePvtx_
 
edm::EDGetTokenT< reco::VertexCollectionvertexProducer_
 
double vtxChi2_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
typedef CacheContexts< T... > CacheTypes
 
typedef CacheTypes::GlobalCache GlobalCache
 
typedef AbilityChecker< T... > HasAbility
 
typedef CacheTypes::LuminosityBlockCache LuminosityBlockCache
 
typedef LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCacheLuminosityBlockContext
 
typedef CacheTypes::LuminosityBlockSummaryCache LuminosityBlockSummaryCache
 
typedef CacheTypes::RunCache RunCache
 
typedef RunContextT< RunCache, GlobalCacheRunContext
 
typedef CacheTypes::RunSummaryCache RunSummaryCache
 

Detailed Description

$Id:

Authors
H. Liu, UC of Riverside US, N. Marinelli Univ of Notre Dame

Description: Produces converted photon objects using default track collections

Implementation: <Notes on="" implementation>="">

Definition at line 64 of file ConversionProducer.h.

Member Typedef Documentation

Definition at line 79 of file ConversionProducer.h.

typedef std::vector<Point> ConversionProducer::PointCollection
private

Definition at line 80 of file ConversionProducer.h.

Constructor & Destructor Documentation

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

Definition at line 67 of file ConversionProducer.cc.

References algoName_, allowD0_, allowDeltaCot_, allowDeltaPhi_, allowMinApproach_, allowOppCharge_, allowSingleLeg_, allowTrackBC_, allowVertex_, bcBarrelCollection_, bcEndcapCollection_, bypassPreselEcal_, bypassPreselEcalEcal_, bypassPreselGsf_, ConvertedPhotonCollection_, d0Cut_, deltaCotTheta_, deltaEta_, deltaPhi_, dEtacutForSCmatching_, dEtaTkBC_, dPhicutForSCmatching_, dPhiTkBC_, dzCut_, energyBC_, energyTotalBC_, edm::ParameterSet::getParameter(), halfWayEta_, maxChi2Left_, maxChi2Right_, maxNumOfTrackInPU_, maxTrackRho_, maxTrackZ_, minApproachHigh_, minApproachLow_, minHitsLeft_, minHitsRight_, minSCEt_, r_cut, rightBC_, scBarrelProducer_, scEndcapProducer_, src_, AlCaHLTBitMon_QueryRunRegistry::string, thettbuilder_, theVertexFinder_, usePvtx_, vertexProducer_, and vtxChi2_.

68  : theVertexFinder_(nullptr)
69 
70 {
71  algoName_ = iConfig.getParameter<std::string>("AlgorithmName");
72 
73  src_ = consumes<edm::View<reco::ConversionTrack> >(iConfig.getParameter<edm::InputTag>("src"));
74 
75  maxNumOfTrackInPU_ = iConfig.getParameter<int>("maxNumOfTrackInPU");
76  maxTrackRho_ = iConfig.getParameter<double>("maxTrackRho");
77  maxTrackZ_ = iConfig.getParameter<double>("maxTrackZ");
78 
79  allowTrackBC_ = iConfig.getParameter<bool>("AllowTrackBC");
80  allowD0_ = iConfig.getParameter<bool>("AllowD0");
81  allowDeltaPhi_ = iConfig.getParameter<bool>("AllowDeltaPhi");
82  allowDeltaCot_ = iConfig.getParameter<bool>("AllowDeltaCot");
83  allowMinApproach_ = iConfig.getParameter<bool>("AllowMinApproach");
84  allowOppCharge_ = iConfig.getParameter<bool>("AllowOppCharge");
85 
86  allowVertex_ = iConfig.getParameter<bool>("AllowVertex");
87 
88  bypassPreselGsf_ = iConfig.getParameter<bool>("bypassPreselGsf");
89  bypassPreselEcal_ = iConfig.getParameter<bool>("bypassPreselEcal");
90  bypassPreselEcalEcal_ = iConfig.getParameter<bool>("bypassPreselEcalEcal");
91 
92  deltaEta_ = iConfig.getParameter<double>("deltaEta");
93 
94  halfWayEta_ = iConfig.getParameter<double>("HalfwayEta"); //open angle to search track matches with BC
95 
96  d0Cut_ = iConfig.getParameter<double>("d0");
97 
98  usePvtx_ = iConfig.getParameter<bool>("UsePvtx"); //if use primary vertices
99 
100  vertexProducer_ = consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("primaryVertexProducer"));
101 
102  //Track-cluster matching eta and phi cuts
103  dEtaTkBC_ = iConfig.getParameter<double>("dEtaTrackBC"); //TODO research on cut endcap/barrel
104  dPhiTkBC_ = iConfig.getParameter<double>("dPhiTrackBC");
105 
107  consumes<edm::View<reco::CaloCluster> >(iConfig.getParameter<edm::InputTag>("bcBarrelCollection"));
109  consumes<edm::View<reco::CaloCluster> >(iConfig.getParameter<edm::InputTag>("bcEndcapCollection"));
110 
111  scBarrelProducer_ = consumes<edm::View<reco::CaloCluster> >(iConfig.getParameter<edm::InputTag>("scBarrelProducer"));
112  scEndcapProducer_ = consumes<edm::View<reco::CaloCluster> >(iConfig.getParameter<edm::InputTag>("scEndcapProducer"));
113 
114  energyBC_ = iConfig.getParameter<double>("EnergyBC"); //BC energy threshold
115  energyTotalBC_ = iConfig.getParameter<double>("EnergyTotalBC"); //BC energy threshold
116  minSCEt_ = iConfig.getParameter<double>("minSCEt"); //super cluster energy threshold
117  dEtacutForSCmatching_ = iConfig.getParameter<double>(
118  "dEtacutForSCmatching"); // dEta between conversion momentum direction and SC position
119  dPhicutForSCmatching_ = iConfig.getParameter<double>(
120  "dPhicutForSCmatching"); // dPhi between conversion momentum direction and SC position
121 
122  //Track cuts on left right track: at least one leg reaches ECAL
123  //Left track: must exist, must reach Ecal and match BC, so loose cut on Chi2 and tight on hits
124  //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
125  maxChi2Left_ = iConfig.getParameter<double>("MaxChi2Left");
126  maxChi2Right_ = iConfig.getParameter<double>("MaxChi2Right");
127  minHitsLeft_ = iConfig.getParameter<int>("MinHitsLeft");
128  minHitsRight_ = iConfig.getParameter<int>("MinHitsRight");
129 
130  //Track Open angle cut on delta cot(theta) and delta phi
131  deltaCotTheta_ = iConfig.getParameter<double>("DeltaCotTheta");
132  deltaPhi_ = iConfig.getParameter<double>("DeltaPhi");
133  minApproachLow_ = iConfig.getParameter<double>("MinApproachLow");
134  minApproachHigh_ = iConfig.getParameter<double>("MinApproachHigh");
135 
136  // if allow single track collection, by default False
137  allowSingleLeg_ = iConfig.getParameter<bool>("AllowSingleLeg");
138  rightBC_ = iConfig.getParameter<bool>("AllowRightBC");
139 
140  //track inner position dz cut, need RECO
141  dzCut_ = iConfig.getParameter<double>("dz");
142  //track analytical cross cut
143  r_cut = iConfig.getParameter<double>("rCut");
144  vtxChi2_ = iConfig.getParameter<double>("vtxChi2");
145 
147 
148  thettbuilder_ = nullptr;
149 
150  //output
151  ConvertedPhotonCollection_ = iConfig.getParameter<std::string>("convertedPhotonCollection");
152 
153  produces<reco::ConversionCollection>(ConvertedPhotonCollection_);
154 }
T getParameter(std::string const &) const
unsigned int maxNumOfTrackInPU_
ConversionVertexFinder * theVertexFinder_
edm::EDGetTokenT< reco::VertexCollection > vertexProducer_
const TransientTrackBuilder * thettbuilder_
edm::EDGetTokenT< edm::View< reco::CaloCluster > > scEndcapProducer_
edm::EDGetTokenT< edm::View< reco::CaloCluster > > bcBarrelCollection_
edm::EDGetTokenT< edm::View< reco::CaloCluster > > bcEndcapCollection_
std::string ConvertedPhotonCollection_
edm::EDGetTokenT< edm::View< reco::ConversionTrack > > src_
edm::EDGetTokenT< edm::View< reco::CaloCluster > > scBarrelProducer_
ConversionProducer::~ConversionProducer ( )
override

Definition at line 156 of file ConversionProducer.cc.

References theVertexFinder_.

156  {
157  // do anything here that needs to be done at desctruction time
158  // (e.g. close files, deallocate resources etc.)
159  delete theVertexFinder_;
160 }
ConversionVertexFinder * theVertexFinder_

Member Function Documentation

void ConversionProducer::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 
)
private

match the track pair with a SC. If at least one track matches, store the SC

Definition at line 274 of file ConversionProducer.cc.

References ecalcalib_dqm_sourceclient-live_cfg::algo, reco::TrackBase::algo(), reco::Conversion::algoByName(), algoName_, allowOppCharge_, allowTrackBC_, muonTagProbeFilters_cff::allTracks, reco::Conversion::arbitratedEcalSeeded, tkConvValidator_cfi::arbitratedEcalSeeded, reco::Conversion::arbitratedMerged, tkConvValidator_cfi::arbitratedMerged, reco::Conversion::arbitratedMergedEcalGeneral, TransientTrackBuilder::build(), bypassPreselEcal_, bypassPreselEcalEcal_, bypassPreselGsf_, edm::RefToBase< T >::castTo(), reco::TrackBase::charge(), checkPhi(), checkVertex(), reco::Vertex::chi2(), ChiSquaredProbability(), deltaEta_, energyTotalBC_, reco::Track::extra(), reco::Conversion::generalTracksOnly, tkConvValidator_cfi::generalTracksOnly, edm::EventSetup::get(), edm::RefToBase< T >::get(), getMatchedBC(), getTrackImpactPosition(), reco::TrackBase::gsf, reco::Conversion::gsfTracksOpenOnly, reco::Conversion::highPurity, reco::Track::innerMomentum(), reco::Track::innerPosition(), reco::TrackBase::inOutEcalSeededConv, edm::Ref< C, T, F >::isNonnull(), reco::Vertex::isValid(), matchingSC(), reco::Vertex::ndof(), ConversionHitChecker::nHitsBeforeVtx(), ConversionHitChecker::nSharedHits(), reco::Track::outerMomentum(), reco::TrackBase::outInEcalSeededConv, preselectTrackPair(), edm::ESHandle< T >::product(), findQualityFiles::rr, reco::Conversion::setMatchingSuperCluster(), reco::Conversion::setQuality(), thettbuilder_, toFConverterP(), toFConverterV(), trackD0Cut(), trackQualityFilter(), and vtxChi2_.

Referenced by produce().

280  {
281  edm::ESHandle<TrackerGeometry> trackerGeomHandle;
282  edm::ESHandle<MagneticField> magFieldHandle;
283  iSetup.get<TrackerDigiGeometryRecord>().get(trackerGeomHandle);
284  iSetup.get<IdealMagneticFieldRecord>().get(magFieldHandle);
285 
286  const TrackerGeometry* trackerGeom = trackerGeomHandle.product();
287  const MagneticField* magField = magFieldHandle.product();
288 
289  // std::vector<math::XYZPointF> trackImpactPosition;
290  // trackImpactPosition.reserve(allTracks.size());//track impact position at ECAL
291  // std::vector<bool> trackValidECAL;//Does this track reach ECAL basic cluster (reach ECAL && match with BC)
292  // trackValidECAL.assign(allTracks.size(), false);
293  //
294  // std::vector<reco::CaloClusterPtr> trackMatchedBC;
295  // reco::CaloClusterPtr empty_bc;
296  // trackMatchedBC.assign(allTracks.size(), empty_bc);//TODO find a better way to avoid copy constructor
297  //
298  // std::vector<int> bcHandleId;//the associated BC handle id, -1 invalid, 0 barrel 1 endcap
299  // bcHandleId.assign(allTracks.size(), -1);
300 
301  // not used std::multimap<double, int> trackInnerEta;//Track innermost state Eta map to TrackRef index, to be used in track pair sorting
302 
303  std::map<edm::Ptr<reco::ConversionTrack>, math::XYZPointF> trackImpactPosition;
304  std::map<edm::Ptr<reco::ConversionTrack>, reco::CaloClusterPtr> trackMatchedBC;
305 
306  ConversionHitChecker hitChecker;
307 
308  //2 propagate all tracks into ECAL, record its eta and phi
309 
310  for (std::multimap<float, edm::Ptr<reco::ConversionTrack> >::const_iterator tk_ref = allTracks.begin();
311  tk_ref != allTracks.end();
312  ++tk_ref) {
313  const reco::Track* tk = tk_ref->second->trackRef().get();
314 
315  //check impact position then match with BC
316  math::XYZPointF ew;
317  if (getTrackImpactPosition(tk, trackerGeom, magField, ew)) {
318  trackImpactPosition[tk_ref->second] = ew;
319 
320  reco::CaloClusterPtr closest_bc; //the closest matching BC to track
321 
322  if (getMatchedBC(basicClusterPtrs, ew, closest_bc)) {
323  trackMatchedBC[tk_ref->second] = closest_bc;
324  }
325  }
326  }
327 
328  //3. pair up tracks:
329  //TODO it is k-Closest pair of point problem
330  //std::cout << " allTracks.size() " << allTracks.size() << std::endl;
331  for (std::multimap<float, edm::Ptr<reco::ConversionTrack> >::const_iterator ll = allTracks.begin();
332  ll != allTracks.end();
333  ++ll) {
334  bool track1HighPurity = true;
335  //std::cout << " Loop on allTracks " << std::endl;
336  const edm::RefToBase<reco::Track>& left = ll->second->trackRef();
337 
338  //TODO: This is a workaround, should be fixed with a proper function in the TTBuilder
339  //(Note that the TrackRef and GsfTrackRef versions of the constructor are needed
340  // to properly get refit tracks in the output vertex)
341  reco::TransientTrack ttk_l;
342  if (dynamic_cast<const reco::GsfTrack*>(left.get())) {
343  ttk_l = thettbuilder_->build(left.castTo<reco::GsfTrackRef>());
344  } else {
345  ttk_l = thettbuilder_->build(left.castTo<reco::TrackRef>());
346  }
347 
349  // if ((allowTrackBC_ && !trackValidECAL[ll-allTracks.begin()]) )//this Left leg should have valid BC
350  // continue;
351 
352  if (the_pvtx.isValid()) {
353  if (!(trackD0Cut(left, the_pvtx)))
354  track1HighPurity = false;
355  } else {
356  if (!(trackD0Cut(left)))
357  track1HighPurity = false;
358  }
359 
360  std::vector<int> right_candidates; //store all right legs passed the cut (theta/approach and ref pair)
361  std::vector<double> right_candidate_theta, right_candidate_approach;
362  std::vector<std::pair<bool, reco::Vertex> > vertex_candidates;
363 
364  //inner loop only over tracks between eta and eta + deltaEta of the first track
365  float etasearch = ll->first + deltaEta_;
366  std::multimap<float, edm::Ptr<reco::ConversionTrack> >::const_iterator rr = ll;
367  ++rr;
368  for (; rr != allTracks.lower_bound(etasearch); ++rr) {
369  bool track2HighPurity = true;
370  bool highPurityPair = true;
371 
372  const edm::RefToBase<reco::Track>& right = rr->second->trackRef();
373 
374  //TODO: This is a workaround, should be fixed with a proper function in the TTBuilder
375  reco::TransientTrack ttk_r;
376  if (dynamic_cast<const reco::GsfTrack*>(right.get())) {
377  ttk_r = thettbuilder_->build(right.castTo<reco::GsfTrackRef>());
378  } else {
379  ttk_r = thettbuilder_->build(right.castTo<reco::TrackRef>());
380  }
381  //std::cout << " This track is " << right->algoName() << std::endl;
382 
383  //all vertexing preselection should go here
384 
385  //check for opposite charge
386  if (allowOppCharge_ && (left->charge() * right->charge() > 0))
387  continue; //same sign, reject pair
388 
390  //if ( (allowTrackBC_ && !trackValidECAL[rr-allTracks.begin()] && rightBC_) )// if right track matches ECAL
391  // continue;
392 
393  double approachDist = -999.;
394  //apply preselection to track pair, overriding preselection for gsf+X or ecalseeded+X pairs if so configured
395  bool preselected = preselectTrackPair(ttk_l, ttk_r, approachDist);
396  preselected = preselected || (bypassPreselGsf_ &&
397  (left->algo() == reco::TrackBase::gsf || right->algo() == reco::TrackBase::gsf));
398  preselected = preselected || (bypassPreselEcal_ && (left->algo() == reco::TrackBase::outInEcalSeededConv ||
402  preselected = preselected || (bypassPreselEcalEcal_ &&
407 
408  if (!preselected) {
409  continue;
410  }
411 
412  //do the actual vertex fit
413  reco::Vertex theConversionVertex; //by default it is invalid
414  bool goodVertex = checkVertex(ttk_l, ttk_r, magField, theConversionVertex);
415 
416  //bail as early as possible in case the fit didn't return a good vertex
417  if (!goodVertex) {
418  continue;
419  }
420 
421  //track pair pass the quality cut
422  if (!((trackQualityFilter(left, true) && trackQualityFilter(right, false)) ||
423  (trackQualityFilter(left, false) && trackQualityFilter(right, true)))) {
424  highPurityPair = false;
425  }
426 
427  if (the_pvtx.isValid()) {
428  if (!(trackD0Cut(right, the_pvtx)))
429  track2HighPurity = false;
430  } else {
431  if (!(trackD0Cut(right)))
432  track2HighPurity = false;
433  }
434 
435  //if all cuts passed, go ahead to make conversion candidates
436  std::vector<edm::RefToBase<reco::Track> > trackPairRef;
437  trackPairRef.push_back(left); //left track
438  trackPairRef.push_back(right); //right track
439 
440  std::vector<math::XYZVectorF> trackPin;
441  std::vector<math::XYZVectorF> trackPout;
442  std::vector<math::XYZPointF> trackInnPos;
443  std::vector<uint8_t> nHitsBeforeVtx;
444  std::vector<Measurement1DFloat> dlClosestHitToVtx;
445 
446  if (left->extra().isNonnull() && right->extra().isNonnull()) { //only available on TrackExtra
447  trackInnPos.push_back(toFConverterP(left->innerPosition()));
448  trackInnPos.push_back(toFConverterP(right->innerPosition()));
449  trackPin.push_back(toFConverterV(left->innerMomentum()));
450  trackPin.push_back(toFConverterV(right->innerMomentum()));
451  trackPout.push_back(toFConverterV(left->outerMomentum()));
452  trackPout.push_back(toFConverterV(right->outerMomentum()));
453  auto leftWrongHits = hitChecker.nHitsBeforeVtx(*left->extra(), theConversionVertex);
454  auto rightWrongHits = hitChecker.nHitsBeforeVtx(*right->extra(), theConversionVertex);
455  nHitsBeforeVtx.push_back(leftWrongHits.first);
456  nHitsBeforeVtx.push_back(rightWrongHits.first);
457  dlClosestHitToVtx.push_back(leftWrongHits.second);
458  dlClosestHitToVtx.push_back(rightWrongHits.second);
459  }
460 
461  uint8_t nSharedHits = hitChecker.nSharedHits(*left.get(), *right.get());
462 
463  //if using kinematic fit, check with chi2 post cut
464  if (theConversionVertex.isValid()) {
465  const float chi2Prob = ChiSquaredProbability(theConversionVertex.chi2(), theConversionVertex.ndof());
466  if (chi2Prob < vtxChi2_)
467  highPurityPair = false;
468  }
469 
470  //std::cout << " highPurityPair after vertex cut " << highPurityPair << std::endl;
471  std::vector<math::XYZPointF> trkPositionAtEcal;
472  std::vector<reco::CaloClusterPtr> matchingBC;
473 
474  if (allowTrackBC_) { //TODO find out the BC ptrs if not doing matching, otherwise, leave it empty
475  //const int lbc_handle = bcHandleId[ll-allTracks.begin()],
476  // rbc_handle = bcHandleId[rr-allTracks.begin()];
477 
478  std::map<edm::Ptr<reco::ConversionTrack>, math::XYZPointF>::const_iterator trackImpactPositionLeft =
479  trackImpactPosition.find(ll->second);
480  std::map<edm::Ptr<reco::ConversionTrack>, math::XYZPointF>::const_iterator trackImpactPositionRight =
481  trackImpactPosition.find(rr->second);
482  std::map<edm::Ptr<reco::ConversionTrack>, reco::CaloClusterPtr>::const_iterator trackMatchedBCLeft =
483  trackMatchedBC.find(ll->second);
484  std::map<edm::Ptr<reco::ConversionTrack>, reco::CaloClusterPtr>::const_iterator trackMatchedBCRight =
485  trackMatchedBC.find(rr->second);
486 
487  if (trackImpactPositionLeft != trackImpactPosition.end()) {
488  trkPositionAtEcal.push_back(trackImpactPositionLeft->second); //left track
489  } else {
490  trkPositionAtEcal.push_back(math::XYZPointF()); //left track
491  }
492  if (trackImpactPositionRight != trackImpactPosition.end()) { //second track ECAL position may be invalid
493  trkPositionAtEcal.push_back(trackImpactPositionRight->second);
494  }
495 
496  double total_e_bc = 0.;
497  if (trackMatchedBCLeft != trackMatchedBC.end()) {
498  matchingBC.push_back(trackMatchedBCLeft->second); //left track
499  total_e_bc += trackMatchedBCLeft->second->energy();
500  } else {
501  matchingBC.push_back(reco::CaloClusterPtr()); //left track
502  }
503  if (trackMatchedBCRight != trackMatchedBC.end()) { //second track ECAL position may be invalid
504  matchingBC.push_back(trackMatchedBCRight->second);
505  total_e_bc += trackMatchedBCRight->second->energy();
506  }
507 
508  if (total_e_bc < energyTotalBC_) {
509  highPurityPair = false;
510  }
511  }
512  //signature cuts, then check if vertex, then post-selection cuts
513  highPurityPair = highPurityPair && track1HighPurity && track2HighPurity && goodVertex &&
514  checkPhi(left, right, trackerGeom, magField, theConversionVertex);
515 
517  /*
518  for ( std::vector<edm::RefToBase<reco::Track> >::iterator iTk=trackPairRef.begin(); iTk!=trackPairRef.end(); iTk++) {
519  math::XYZPointF impPos;
520  if ( getTrackImpactPosition(*iTk, trackerGeom, magField, impPos) ) {
521 
522  }
523 
524  }
525  */
526 
527  const float minAppDist = approachDist;
529  float dummy = 0;
531  reco::Conversion newCandidate(scPtrVec,
532  trackPairRef,
533  trkPositionAtEcal,
534  theConversionVertex,
535  matchingBC,
536  minAppDist,
537  trackInnPos,
538  trackPin,
539  trackPout,
540  nHitsBeforeVtx,
541  dlClosestHitToVtx,
542  nSharedHits,
543  dummy,
544  algo);
545  // Fill in scPtrVec with the macthing SC
546  if (matchingSC(superClusterPtrs, newCandidate, scPtrVec))
547  newCandidate.setMatchingSuperCluster(scPtrVec);
548 
549  newCandidate.setQuality(reco::Conversion::highPurity, highPurityPair);
550  bool generalTracksOnly = ll->second->isTrackerOnly() && rr->second->isTrackerOnly() &&
551  !dynamic_cast<const reco::GsfTrack*>(ll->second->trackRef().get()) &&
552  !dynamic_cast<const reco::GsfTrack*>(rr->second->trackRef().get());
553  bool gsfTracksOpenOnly = ll->second->isGsfTrackOpen() && rr->second->isGsfTrackOpen();
554  bool arbitratedEcalSeeded = ll->second->isArbitratedEcalSeeded() && rr->second->isArbitratedEcalSeeded();
555  bool arbitratedMerged = ll->second->isArbitratedMerged() && rr->second->isArbitratedMerged();
556  bool arbitratedMergedEcalGeneral =
557  ll->second->isArbitratedMergedEcalGeneral() && rr->second->isArbitratedMergedEcalGeneral();
558 
559  newCandidate.setQuality(reco::Conversion::generalTracksOnly, generalTracksOnly);
560  newCandidate.setQuality(reco::Conversion::gsfTracksOpenOnly, gsfTracksOpenOnly);
561  newCandidate.setQuality(reco::Conversion::arbitratedEcalSeeded, arbitratedEcalSeeded);
562  newCandidate.setQuality(reco::Conversion::arbitratedMerged, arbitratedMerged);
563  newCandidate.setQuality(reco::Conversion::arbitratedMergedEcalGeneral, arbitratedMergedEcalGeneral);
564 
565  outputConvPhotonCollection.push_back(newCandidate);
566  }
567  }
568 }
value_type const * get() const
Definition: RefToBase.h:209
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
bool getTrackImpactPosition(const reco::Track *tk_ref, const TrackerGeometry *trackerGeom, const MagneticField *magField, math::XYZPointF &ew)
const TrackExtraRef & extra() const
reference to "extra" object
Definition: Track.h:139
const TransientTrackBuilder * thettbuilder_
bool isValid() const
Tells whether the vertex is valid.
Definition: Vertex.h:71
reco::TransientTrack build(const reco::Track *p) const
static ConversionAlgorithm algoByName(const std::string &name)
Definition: Conversion.cc:139
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float > > XYZPointF
point in space with cartesian internal representation
Definition: Point3D.h:10
std::pair< uint8_t, Measurement1DFloat > nHitsBeforeVtx(const reco::TrackExtra &track, const reco::Vertex &vtx, float sigmaTolerance=3.0) const
bool matchingSC(const std::multimap< double, reco::CaloClusterPtr > &scMap, reco::Conversion &conv, reco::CaloClusterPtrVector &mSC)
const math::XYZPoint & innerPosition() const
position of the innermost hit
Definition: Track.h:56
TrackAlgorithm algo() const
Definition: TrackBase.h:526
math::XYZVectorF toFConverterV(const math::XYZVector &val)
double chi2() const
chi-squares
Definition: Vertex.h:102
float ChiSquaredProbability(double chiSquared, double nrDOF)
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:109
bool trackQualityFilter(const edm::RefToBase< reco::Track > &ref, bool isLeft)
const math::XYZVector & outerMomentum() const
momentum vector at the outermost hit position
Definition: Track.h:65
bool trackD0Cut(const edm::RefToBase< reco::Track > &ref)
REF castTo() const
Definition: RefToBase.h:257
const math::XYZVector & innerMomentum() const
momentum vector at the innermost hit position
Definition: Track.h:59
T get() const
Definition: EventSetup.h:73
int charge() const
track electric charge
Definition: TrackBase.h:575
math::XYZPointF toFConverterP(const math::XYZPoint &val)
bool preselectTrackPair(const reco::TransientTrack &ttk_l, const reco::TransientTrack &ttk_r, double &appDist)
T const * product() const
Definition: ESHandle.h:86
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)
bool checkVertex(const reco::TransientTrack &ttk_l, const reco::TransientTrack &ttk_r, const MagneticField *magField, reco::Vertex &the_vertex)
void ConversionProducer::buildSuperAndBasicClusterGeoMap ( const edm::Event iEvent,
std::multimap< double, reco::CaloClusterPtr > &  basicClusterPtrs,
std::multimap< double, reco::CaloClusterPtr > &  superClusterPtrs 
)
private

Definition at line 225 of file ConversionProducer.cc.

References bcBarrelCollection_, bcEndcapCollection_, SimCluster::energy(), energyBC_, SimCluster::eta(), edm::Event::getByToken(), patZpeak::handle, edm::HandleBase::isValid(), minSCEt_, SimDataFormats::CaloAnalysis::sc, scBarrelProducer_, and scEndcapProducer_.

Referenced by produce().

227  {
228  // Get the Super Cluster collection in the Barrel
230  iEvent.getByToken(scBarrelProducer_, scBarrelHandle);
231  if (!scBarrelHandle.isValid()) {
232  edm::LogError("ConvertedPhotonProducer") << "Error! Can't get the barrel superclusters!";
233  }
234 
235  // Get the Super Cluster collection in the Endcap
237  iEvent.getByToken(scEndcapProducer_, scEndcapHandle);
238  if (!scEndcapHandle.isValid()) {
239  edm::LogError("ConvertedPhotonProducer") << "Error! Can't get the endcap superclusters!";
240  }
241 
243  edm::Handle<edm::View<reco::CaloCluster> > bcEndcapHandle; //TODO check cluster type if BasicCluster or PFCluster
244 
245  iEvent.getByToken(bcBarrelCollection_, bcBarrelHandle);
246  if (!bcBarrelHandle.isValid()) {
247  edm::LogError("ConvertedPhotonProducer") << "Error! Can't get the barrel basic clusters!";
248  }
249 
250  iEvent.getByToken(bcEndcapCollection_, bcEndcapHandle);
251  if (!bcEndcapHandle.isValid()) {
252  edm::LogError("ConvertedPhotonProducer") << "Error! Can't get the endcap basic clusters!";
253  }
254 
255  if (bcBarrelHandle.isValid()) {
256  for (auto const& handle : {bcBarrelHandle, bcEndcapHandle}) {
257  for (auto const& bc : handle->ptrs()) {
258  if (bc->energy() > energyBC_)
259  basicClusterPtrs.emplace(bc->position().eta(), bc);
260  }
261  }
262  }
263 
264  if (scBarrelHandle.isValid()) {
265  for (auto const& handle : {scBarrelHandle, scEndcapHandle}) {
266  for (auto const& sc : handle->ptrs()) {
267  if (sc->energy() > minSCEt_)
268  superClusterPtrs.emplace(sc->position().eta(), sc);
269  }
270  }
271  }
272 }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
float energy() const
Energy. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:104
edm::EDGetTokenT< edm::View< reco::CaloCluster > > scEndcapProducer_
edm::EDGetTokenT< edm::View< reco::CaloCluster > > bcBarrelCollection_
bool isValid() const
Definition: HandleBase.h:70
edm::EDGetTokenT< edm::View< reco::CaloCluster > > bcEndcapCollection_
edm::EDGetTokenT< edm::View< reco::CaloCluster > > scBarrelProducer_
float eta() const
Momentum pseudorapidity. Note this is taken from the simtrack before the calorimeter.
Definition: SimCluster.h:148
bool ConversionProducer::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 
)
private

Definition at line 700 of file ConversionProducer.cc.

References allowDeltaPhi_, anyDirection, reco::deltaPhi(), deltaPhi_, HLT_2018_cff::dPhi, reco::Track::extra(), TrajectoryStateOnSurface::globalDirection(), reco::Track::innerMomentum(), trajectoryStateTransform::innerStateOnSurface(), edm::Ref< C, T, F >::isNonnull(), TrajectoryStateOnSurface::isValid(), reco::Vertex::isValid(), PV3DBase< T, PVType, FrameType >::phi(), reco::Vertex::position(), Propagator::propagate(), and makeMuonMisalignmentScenario::rot.

Referenced by buildCollection().

704  {
705  if (!allowDeltaPhi_)
706  return true;
707  //if track has innermost momentum, check with innermost phi
708  //if track also has valid vertex, propagate to vertex then calculate phi there
709  //if track has no innermost momentum, just return true, because track->phi() makes no sense
710  if (tk_l->extra().isNonnull() && tk_r->extra().isNonnull()) {
711  double iphi1 = tk_l->innerMomentum().phi(), iphi2 = tk_r->innerMomentum().phi();
712  if (vtx.isValid()) {
713  PropagatorWithMaterial propag(anyDirection, 0.000511, magField);
714 
715  double recoPhoR = vtx.position().Rho();
718  recoPhoR,
719  Surface::PositionType(0, 0, 0),
720  rot,
722  recoPhoR - 0.001, recoPhoR + 0.001, -fabs(vtx.position().z()), fabs(vtx.position().z()))));
724  Surface::PositionType(0, 0, vtx.position().z()), rot, new SimpleDiskBounds(0, recoPhoR, -0.001, 0.001)));
725 
726  const TrajectoryStateOnSurface myTSOS1 =
727  trajectoryStateTransform::innerStateOnSurface(*tk_l, *trackerGeom, magField);
728  const TrajectoryStateOnSurface myTSOS2 =
729  trajectoryStateTransform::innerStateOnSurface(*tk_r, *trackerGeom, magField);
730  TrajectoryStateOnSurface stateAtVtx1, stateAtVtx2;
731  stateAtVtx1 = propag.propagate(myTSOS1, *theBarrel_);
732  if (!stateAtVtx1.isValid()) {
733  stateAtVtx1 = propag.propagate(myTSOS1, *theDisk_);
734  }
735  if (stateAtVtx1.isValid()) {
736  iphi1 = stateAtVtx1.globalDirection().phi();
737  }
738  stateAtVtx2 = propag.propagate(myTSOS2, *theBarrel_);
739  if (!stateAtVtx2.isValid()) {
740  stateAtVtx2 = propag.propagate(myTSOS2, *theDisk_);
741  }
742  if (stateAtVtx2.isValid()) {
743  iphi2 = stateAtVtx2.globalDirection().phi();
744  }
745  }
746  const double dPhi = reco::deltaPhi(iphi1, iphi2);
747  return (fabs(dPhi) < deltaPhi_);
748  } else {
749  return true;
750  }
751 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
const TrackExtraRef & extra() const
reference to "extra" object
Definition: Track.h:139
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
Cylinder BoundCylinder
Definition: BoundCylinder.h:17
Disk BoundDisk
Definition: BoundDisk.h:54
const math::XYZVector & innerMomentum() const
momentum vector at the innermost hit position
Definition: Track.h:59
GlobalVector globalDirection() const
TrajectoryStateOnSurface innerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field, bool withErr=true)
bool ConversionProducer::checkTrackPair ( const std::pair< edm::RefToBase< reco::Track >, reco::CaloClusterPtr > &  ll,
const std::pair< edm::RefToBase< reco::Track >, reco::CaloClusterPtr > &  rr 
)
private

Definition at line 808 of file ConversionProducer.cc.

References allowTrackBC_, energyTotalBC_, edm::Ptr< T >::isNonnull(), rightBC_, and findQualityFiles::rr.

809  {
810  const reco::CaloClusterPtr& bc_l = ll.second; //can be null, so check isNonnull()
811  const reco::CaloClusterPtr& bc_r = rr.second;
812 
813  //The cuts should be ordered by considering if takes time and if cuts off many fakes
814  if (allowTrackBC_) {
815  //check energy of BC
816  double total_e_bc = 0;
817  if (bc_l.isNonnull())
818  total_e_bc += bc_l->energy();
819  if (rightBC_)
820  if (bc_r.isNonnull())
821  total_e_bc += bc_r->energy();
822 
823  if (total_e_bc < energyTotalBC_)
824  return false;
825  }
826 
827  return true;
828 }
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:146
bool ConversionProducer::checkVertex ( const reco::TransientTrack ttk_l,
const reco::TransientTrack ttk_r,
const MagneticField magField,
reco::Vertex the_vertex 
)
private

Definition at line 831 of file ConversionProducer.cc.

References newFWLiteAna::found, ConversionVertexFinder::run(), and theVertexFinder_.

Referenced by buildCollection().

834  {
835  bool found = false;
836 
837  std::vector<reco::TransientTrack> pair;
838  pair.push_back(ttk_l);
839  pair.push_back(ttk_r);
840 
841  found = theVertexFinder_->run(pair, the_vertex);
842 
843  return found;
844 }
ConversionVertexFinder * theVertexFinder_
TransientVertex run(const std::vector< reco::TransientTrack > &pair)
double ConversionProducer::etaTransformation ( float  EtaParticle,
float  Zvertex 
)
private

Definition at line 846 of file ConversionProducer.cc.

References ETA, etaBarrelEndcap, dqm-mbProfile::log, PI, R_ECAL, funct::tan(), and Z_Endcap.

Referenced by matchingSC().

846  {
847  //---Definitions
848  const float PI = 3.1415927;
849 
850  //---Definitions for ECAL
851  const float R_ECAL = 136.5;
852  const float Z_Endcap = 328.0;
853  const float etaBarrelEndcap = 1.479;
854 
855  //---ETA correction
856 
857  float Theta = 0.0;
858  float ZEcal = R_ECAL * sinh(EtaParticle) + Zvertex;
859 
860  if (ZEcal != 0.0)
861  Theta = atan(R_ECAL / ZEcal);
862  if (Theta < 0.0)
863  Theta = Theta + PI;
864  double ETA = -log(tan(0.5 * Theta));
865 
866  if (fabs(ETA) > etaBarrelEndcap) {
867  float Zend = Z_Endcap;
868  if (EtaParticle < 0.0)
869  Zend = -Zend;
870  float Zlen = Zend - Zvertex;
871  float RR = Zlen / sinh(EtaParticle);
872  Theta = atan(RR / Zend);
873  if (Theta < 0.0)
874  Theta = Theta + PI;
875  ETA = -log(tan(0.5 * Theta));
876  }
877  //---Return the result
878  return ETA;
879  //---end
880 }
#define ETA
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
static float etaBarrelEndcap
#define PI
Definition: QcdUeDQM.h:37
static float Z_Endcap
static float R_ECAL
bool ConversionProducer::getMatchedBC ( const std::multimap< double, reco::CaloClusterPtr > &  bcMap,
const math::XYZPointF trackImpactPosition,
reco::CaloClusterPtr closestBC 
)
private

Definition at line 668 of file ConversionProducer.cc.

References HLT_2018_cff::delta_eta, HLT_2018_cff::delta_phi, reco::deltaPhi(), dEtaTkBC_, dPhiTkBC_, halfWayEta_, and HLT_2018_cff::min_eta.

Referenced by buildCollection().

670  {
671  const double track_eta = trackImpactPosition.eta();
672  const double track_phi = trackImpactPosition.phi();
673 
674  double min_eta = 999., min_phi = 999.;
675  reco::CaloClusterPtr closest_bc;
676  for (std::multimap<double, reco::CaloClusterPtr>::const_iterator bc = bcMap.lower_bound(track_eta - halfWayEta_);
677  bc != bcMap.upper_bound(track_eta + halfWayEta_);
678  ++bc) { //use eta map to select possible BC collection then loop in
679  const reco::CaloClusterPtr& ebc = bc->second;
680  const double delta_eta = track_eta - (ebc->position().eta());
681  const double delta_phi = reco::deltaPhi(track_phi, (ebc->position().phi()));
682  if (fabs(delta_eta) < dEtaTkBC_ && fabs(delta_phi) < dPhiTkBC_) {
683  if (fabs(min_eta) > fabs(delta_eta) && fabs(min_phi) > fabs(delta_phi)) { //take the closest to track BC
684  min_eta = delta_eta;
685  min_phi = delta_phi;
686  closest_bc = bc->second;
687  //TODO check if min_eta>delta_eta but min_phi<delta_phi
688  }
689  }
690  }
691 
692  if (min_eta < 999.) {
693  closestBC = closest_bc;
694  return true;
695  } else
696  return false;
697 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
bool ConversionProducer::getTrackImpactPosition ( const reco::Track tk_ref,
const TrackerGeometry trackerGeom,
const MagneticField magField,
math::XYZPointF ew 
)
private

Definition at line 595 of file ConversionProducer.cc.

References alongMomentum, barrelHalfLength, barrelRadius, endcapRadius, endcapZ, geometryDiff::epsilon, PV3DBase< T, PVType, FrameType >::eta(), f, TrajectoryStateOnSurface::globalPosition(), TrajectoryStateOnSurface::isValid(), trajectoryStateTransform::outerStateOnSurface(), Propagator::propagate(), makeMuonMisalignmentScenario::rot, and PV3DBase< T, PVType, FrameType >::z().

Referenced by buildCollection().

598  {
599  PropagatorWithMaterial propag(alongMomentum, 0.000511, magField);
600 
602  129.f, GlobalPoint(0., 0., 0.), TkRotation<float>(), new SimpleCylinderBounds(129, 129, -320.5, 320.5)));
603  const float epsilon = 0.001;
604  Surface::RotationType rot; // unit rotation matrix
605  const float barrelRadius = 129.f;
606  const float barrelHalfLength = 270.9f;
607  const float endcapRadius = 171.1f;
608  const float endcapZ = 320.5f;
610  barrelRadius,
611  Surface::PositionType(0, 0, 0),
612  rot,
613  new SimpleCylinderBounds(barrelRadius - epsilon, barrelRadius + epsilon, -barrelHalfLength, barrelHalfLength)));
614  ReferenceCountingPointer<BoundDisk> theNegativeEtaEndcap_(new BoundDisk(
615  Surface::PositionType(0, 0, -endcapZ), rot, new SimpleDiskBounds(0, endcapRadius, -epsilon, epsilon)));
616  ReferenceCountingPointer<BoundDisk> thePositiveEtaEndcap_(new BoundDisk(
617  Surface::PositionType(0, 0, endcapZ), rot, new SimpleDiskBounds(0, endcapRadius, -epsilon, epsilon)));
618 
619  //const TrajectoryStateOnSurface myTSOS = trajectoryStateTransform::innerStateOnSurface(*(*ref), *trackerGeom, magField);
620  const TrajectoryStateOnSurface myTSOS =
621  trajectoryStateTransform::outerStateOnSurface(*tk_ref, *trackerGeom, magField);
622  TrajectoryStateOnSurface stateAtECAL;
623  stateAtECAL = propag.propagate(myTSOS, *theBarrel_);
624  if (!stateAtECAL.isValid() || (stateAtECAL.isValid() && fabs(stateAtECAL.globalPosition().eta()) > 1.479f)) {
625  //endcap propagator
626  if (myTSOS.globalPosition().z() > 0.) {
627  stateAtECAL = propag.propagate(myTSOS, *thePositiveEtaEndcap_);
628  } else {
629  stateAtECAL = propag.propagate(myTSOS, *theNegativeEtaEndcap_);
630  }
631  }
632  if (stateAtECAL.isValid()) {
633  ew = stateAtECAL.globalPosition();
634  return true;
635  } else
636  return false;
637 }
TrajectoryStateOnSurface outerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field, bool withErr=true)
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
Cylinder BoundCylinder
Definition: BoundCylinder.h:17
GlobalPoint globalPosition() const
T z() const
Definition: PV3DBase.h:61
double f[11][100]
Disk BoundDisk
Definition: BoundDisk.h:54
T eta() const
Definition: PV3DBase.h:73
bool ConversionProducer::matchingSC ( const std::multimap< double, reco::CaloClusterPtr > &  scMap,
reco::Conversion conv,
reco::CaloClusterPtrVector mSC 
)
private

Definition at line 639 of file ConversionProducer.cc.

References HLT_2018_cff::delta_eta, HLT_2018_cff::delta_phi, reco::deltaPhi(), dEtacutForSCmatching_, dPhicutForSCmatching_, TrackSplittingMonitor_cfi::dphiMin, etaTransformation(), match(), edm::PtrVector< T >::push_back(), reco::Conversion::refittedPairMomentum(), SimDataFormats::CaloAnalysis::sc, and reco::Conversion::zOfPrimaryVertexFromTracks().

Referenced by buildCollection().

642  {
643  // double dRMin=999.;
644  double detaMin = 999.;
645  double dphiMin = 999.;
647  for (std::multimap<double, reco::CaloClusterPtr>::const_iterator scItr = scMap.begin(); scItr != scMap.end();
648  scItr++) {
649  const reco::CaloClusterPtr& sc = scItr->second;
650  const double delta_phi = reco::deltaPhi(aConv.refittedPairMomentum().phi(), sc->phi());
651  double sceta = sc->eta();
652  double conveta = etaTransformation(aConv.refittedPairMomentum().eta(), aConv.zOfPrimaryVertexFromTracks());
653  const double delta_eta = fabs(conveta - sceta);
654  if (fabs(delta_eta) < fabs(detaMin) && fabs(delta_phi) < fabs(dphiMin)) {
655  detaMin = fabs(delta_eta);
656  dphiMin = fabs(delta_phi);
657  match = sc;
658  }
659  }
660 
661  if (fabs(detaMin) < dEtacutForSCmatching_ && fabs(dphiMin) < dPhicutForSCmatching_) {
662  mSC.push_back(match);
663  return true;
664  } else
665  return false;
666 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:149
double etaTransformation(float EtaParticle, float Zvertex)
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
bool ConversionProducer::preselectTrackPair ( const reco::TransientTrack ttk_l,
const reco::TransientTrack ttk_r,
double &  appDist 
)
private

Definition at line 753 of file ConversionProducer.cc.

References funct::abs(), allowDeltaCot_, allowMinApproach_, TangentApproachInRPhi::calculate(), ClosestApproachInRPhi::calculate(), TangentApproachInRPhi::crossingPoint(), ClosestApproachInRPhi::crossingPoint(), deltaCotTheta_, dzCut_, reco::Track::innerMomentum(), reco::TransientTrack::innermostMeasurementState(), maxTrackRho_, maxTrackZ_, minApproachHigh_, PV3DBase< T, PVType, FrameType >::perp(), TangentApproachInRPhi::perpdist(), r_cut, rho, TangentApproachInRPhi::status(), ClosestApproachInRPhi::status(), funct::tan(), reco::TransientTrack::track(), TangentApproachInRPhi::trajectoryParameters(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by buildCollection().

755  {
756  double dCotTheta = 1. / tan(ttk_l.track().innerMomentum().theta()) - 1. / tan(ttk_r.track().innerMomentum().theta());
757  if (allowDeltaCot_ && (std::abs(dCotTheta) > deltaCotTheta_)) {
758  return false;
759  }
760 
761  //non-conversion hypothesis, reject prompt track pairs
762  ClosestApproachInRPhi closest;
764  if (!closest.status()) {
765  return false;
766  }
767 
768  if (closest.crossingPoint().perp() < r_cut) {
769  return false;
770  }
771 
772  //compute tangent point btw tracks (conversion hypothesis)
773  TangentApproachInRPhi tangent;
775  if (!tangent.status()) {
776  return false;
777  }
778 
779  GlobalPoint tangentPoint = tangent.crossingPoint();
780  double rho = tangentPoint.perp();
781 
782  //reject candidates well outside of tracker bounds
783  if (rho > maxTrackRho_) {
784  return false;
785  }
786 
787  if (std::abs(tangentPoint.z()) > maxTrackZ_) {
788  return false;
789  }
790 
791  std::pair<GlobalTrajectoryParameters, GlobalTrajectoryParameters> trajs = tangent.trajectoryParameters();
792 
793  //very large separation in z, no hope
794  if (std::abs(trajs.first.position().z() - trajs.second.position().z()) > dzCut_) {
795  return false;
796  }
797 
798  float minApproach = tangent.perpdist();
799  appDist = minApproach;
800 
801  if (allowMinApproach_ && (minApproach < minApproachLow_ || minApproach > minApproachHigh_)) {
802  return false;
803  }
804 
805  return true;
806 }
T perp() const
Definition: PV3DBase.h:69
std::pair< GlobalTrajectoryParameters, GlobalTrajectoryParameters > trajectoryParameters() const
TrajectoryStateOnSurface innermostMeasurementState() const
bool status() const override
bool calculate(const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb) override
T z() const
Definition: PV3DBase.h:61
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool calculate(const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb) override
const Track & track() const
GlobalPoint crossingPoint() const override
const math::XYZVector & innerMomentum() const
momentum vector at the innermost hit position
Definition: Track.h:59
bool status() const override
GlobalPoint crossingPoint() const override
void ConversionProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 163 of file ConversionProducer.cc.

References buildCollection(), buildSuperAndBasicClusterGeoMap(), ConvertedPhotonCollection_, edm::EventSetup::get(), edm::Event::getByToken(), edm::HandleBase::isValid(), maxNumOfTrackInPU_, eostools::move(), edm::Handle< T >::product(), edm::ESHandle< T >::product(), edm::Event::put(), src_, OrderedSet::t, thettbuilder_, usePvtx_, spclusmultinvestigator_cfi::vertexCollection, and vertexProducer_.

163  {
164  using namespace edm;
165 
166  reco::ConversionCollection outputConvPhotonCollection;
167  auto outputConvPhotonCollection_p = std::make_unique<reco::ConversionCollection>();
168 
169  //std::cout << " ConversionProducer::produce " << std::endl;
170  //Read multiple track input collections
171 
172  edm::Handle<edm::View<reco::ConversionTrack> > trackCollectionHandle;
173  iEvent.getByToken(src_, trackCollectionHandle);
174 
175  //build map of ConversionTracks ordered in eta
176  std::multimap<float, edm::Ptr<reco::ConversionTrack> > convTrackMap;
177  for (auto const& t : trackCollectionHandle->ptrs())
178  convTrackMap.emplace(t->track()->eta(), t);
179 
182  if (usePvtx_) {
183  iEvent.getByToken(vertexProducer_, vertexHandle);
184  if (!vertexHandle.isValid()) {
185  edm::LogError("ConversionProducer") << "Error! Can't get the product primary Vertex Collection "
186  << "\n";
187  usePvtx_ = false;
188  }
189  if (usePvtx_)
190  vertexCollection = *(vertexHandle.product());
191  }
192 
193  edm::ESHandle<TransientTrackBuilder> hTransientTrackBuilder;
194  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder", hTransientTrackBuilder);
195  thettbuilder_ = hTransientTrackBuilder.product();
196 
197  reco::Vertex the_pvtx;
198  //because the priamry vertex is sorted by quality, the first one is the best
199  if (!vertexCollection.empty())
200  the_pvtx = *(vertexCollection.begin());
201 
202  if (trackCollectionHandle->size() > maxNumOfTrackInPU_) {
203  iEvent.put(std::move(outputConvPhotonCollection_p), ConvertedPhotonCollection_);
204  return;
205  }
206 
207  // build Super and Basic cluster geometry map to search in eta bounds for clusters
208  std::multimap<double, reco::CaloClusterPtr> basicClusterPtrs;
209  std::multimap<double, reco::CaloClusterPtr> superClusterPtrs;
210 
211  buildSuperAndBasicClusterGeoMap(iEvent, basicClusterPtrs, superClusterPtrs);
212 
213  buildCollection(iEvent,
214  iSetup,
215  convTrackMap,
216  superClusterPtrs,
217  basicClusterPtrs,
218  the_pvtx,
219  outputConvPhotonCollection); //allow empty basicClusterPtrs
220 
221  outputConvPhotonCollection_p->assign(outputConvPhotonCollection.begin(), outputConvPhotonCollection.end());
222  iEvent.put(std::move(outputConvPhotonCollection_p), ConvertedPhotonCollection_);
223 }
unsigned int maxNumOfTrackInPU_
edm::EDGetTokenT< reco::VertexCollection > vertexProducer_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
const TransientTrackBuilder * thettbuilder_
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)
void buildSuperAndBasicClusterGeoMap(const edm::Event &, std::multimap< double, reco::CaloClusterPtr > &basicClusterPtrs, std::multimap< double, reco::CaloClusterPtr > &superClusterPtrs)
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
std::vector< Conversion > ConversionCollection
collectin of Conversion objects
Definition: ConversionFwd.h:9
bool isValid() const
Definition: HandleBase.h:70
T const * product() const
Definition: Handle.h:69
std::string ConvertedPhotonCollection_
edm::EDGetTokenT< edm::View< reco::ConversionTrack > > src_
HLT enums.
T get() const
Definition: EventSetup.h:73
T const * product() const
Definition: ESHandle.h:86
def move(src, dest)
Definition: eostools.py:511
math::XYZPointF ConversionProducer::toFConverterP ( const math::XYZPoint val)
inlineprivate

Definition at line 179 of file ConversionProducer.h.

Referenced by buildCollection().

179 { return math::XYZPointF(val.x(), val.y(), val.z()); }
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float > > XYZPointF
point in space with cartesian internal representation
Definition: Point3D.h:10
math::XYZVectorF ConversionProducer::toFConverterV ( const math::XYZVector val)
inlineprivate

Definition at line 181 of file ConversionProducer.h.

Referenced by buildCollection().

181 { return math::XYZVectorF(val.x(), val.y(), val.z()); }
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float > > XYZVectorF
spatial vector with cartesian internal representation
Definition: Vector3D.h:16
bool ConversionProducer::trackD0Cut ( const edm::RefToBase< reco::Track > &  ref)
inlineprivate

Definition at line 585 of file ConversionProducer.cc.

References allowD0_, reco::TrackBase::charge(), reco::TrackBase::d0(), d0Cut_, and reco::TrackBase::d0Error().

Referenced by buildCollection().

585  {
586  //NOTE if not allow d0 cut, always true
587  return ((!allowD0_) || !(ref->d0() * ref->charge() / ref->d0Error() < d0Cut_));
588 }
double d0Error() const
error on d0
Definition: TrackBase.h:719
double d0() const
dxy parameter in perigee convention (d0 = -dxy)
Definition: TrackBase.h:590
int charge() const
track electric charge
Definition: TrackBase.h:575
bool ConversionProducer::trackD0Cut ( const edm::RefToBase< reco::Track > &  ref,
const reco::Vertex the_pvtx 
)
inlineprivate

Definition at line 590 of file ConversionProducer.cc.

References allowD0_, reco::TrackBase::charge(), d0Cut_, reco::TrackBase::dxy(), reco::TrackBase::dxyError(), and reco::Vertex::position().

590  {
591  //
592  return ((!allowD0_) || !(-ref->dxy(the_pvtx.position()) * ref->charge() / ref->dxyError() < d0Cut_));
593 }
double dxyError() const
error on dxy
Definition: TrackBase.h:716
const Point & position() const
position
Definition: Vertex.h:113
int charge() const
track electric charge
Definition: TrackBase.h:575
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:587
bool ConversionProducer::trackQualityFilter ( const edm::RefToBase< reco::Track > &  ref,
bool  isLeft 
)
inlineprivate

Definition at line 574 of file ConversionProducer.cc.

References reco::Track::found(), maxChi2Left_, maxChi2Right_, minHitsLeft_, minHitsRight_, and reco::TrackBase::normalizedChi2().

Referenced by buildCollection().

574  {
575  bool pass = true;
576  if (isLeft) {
577  pass = (ref->normalizedChi2() < maxChi2Left_ && ref->found() >= minHitsLeft_);
578  } else {
579  pass = (ref->normalizedChi2() < maxChi2Right_ && ref->found() >= minHitsRight_);
580  }
581 
582  return pass;
583 }
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:572
unsigned short found() const
Number of valid hits on track.
Definition: Track.h:142

Member Data Documentation

std::string ConversionProducer::algoName_
private

Definition at line 77 of file ConversionProducer.h.

Referenced by buildCollection(), and ConversionProducer().

bool ConversionProducer::allowD0_
private

Definition at line 90 of file ConversionProducer.h.

Referenced by ConversionProducer(), and trackD0Cut().

bool ConversionProducer::allowDeltaCot_
private

Definition at line 90 of file ConversionProducer.h.

Referenced by ConversionProducer(), and preselectTrackPair().

bool ConversionProducer::allowDeltaPhi_
private

Definition at line 90 of file ConversionProducer.h.

Referenced by checkPhi(), and ConversionProducer().

bool ConversionProducer::allowMinApproach_
private

Definition at line 90 of file ConversionProducer.h.

Referenced by ConversionProducer(), and preselectTrackPair().

bool ConversionProducer::allowOppCharge_
private

Definition at line 90 of file ConversionProducer.h.

Referenced by buildCollection(), and ConversionProducer().

bool ConversionProducer::allowSingleLeg_
private

Definition at line 124 of file ConversionProducer.h.

Referenced by ConversionProducer().

bool ConversionProducer::allowTrackBC_
private

Definition at line 90 of file ConversionProducer.h.

Referenced by buildCollection(), checkTrackPair(), and ConversionProducer().

bool ConversionProducer::allowVertex_
private

Definition at line 90 of file ConversionProducer.h.

Referenced by ConversionProducer().

edm::EDGetTokenT<edm::View<reco::CaloCluster> > ConversionProducer::bcBarrelCollection_
private

Definition at line 86 of file ConversionProducer.h.

Referenced by buildSuperAndBasicClusterGeoMap(), and ConversionProducer().

edm::EDGetTokenT<edm::View<reco::CaloCluster> > ConversionProducer::bcEndcapCollection_
private

Definition at line 87 of file ConversionProducer.h.

Referenced by buildSuperAndBasicClusterGeoMap(), and ConversionProducer().

bool ConversionProducer::bypassPreselEcal_
private

Definition at line 92 of file ConversionProducer.h.

Referenced by buildCollection(), and ConversionProducer().

bool ConversionProducer::bypassPreselEcalEcal_
private

Definition at line 92 of file ConversionProducer.h.

Referenced by buildCollection(), and ConversionProducer().

bool ConversionProducer::bypassPreselGsf_
private

Definition at line 92 of file ConversionProducer.h.

Referenced by buildCollection(), and ConversionProducer().

std::string ConversionProducer::ConvertedPhotonCollection_
private

Definition at line 88 of file ConversionProducer.h.

Referenced by ConversionProducer(), and produce().

double ConversionProducer::d0Cut_
private

Definition at line 111 of file ConversionProducer.h.

Referenced by ConversionProducer(), and trackD0Cut().

double ConversionProducer::deltaCotTheta_
private

Definition at line 118 of file ConversionProducer.h.

Referenced by ConversionProducer(), and preselectTrackPair().

double ConversionProducer::deltaEta_
private

Definition at line 100 of file ConversionProducer.h.

Referenced by buildCollection(), and ConversionProducer().

double ConversionProducer::deltaPhi_
private

Definition at line 118 of file ConversionProducer.h.

Referenced by checkPhi(), and ConversionProducer().

double ConversionProducer::dEtacutForSCmatching_
private

Definition at line 107 of file ConversionProducer.h.

Referenced by ConversionProducer(), and matchingSC().

double ConversionProducer::dEtaTkBC_
private

Definition at line 113 of file ConversionProducer.h.

Referenced by ConversionProducer(), and getMatchedBC().

double ConversionProducer::dPhicutForSCmatching_
private

Definition at line 108 of file ConversionProducer.h.

Referenced by ConversionProducer(), and matchingSC().

double ConversionProducer::dPhiTkBC_
private

Definition at line 113 of file ConversionProducer.h.

Referenced by ConversionProducer(), and getMatchedBC().

double ConversionProducer::dzCut_
private

Definition at line 112 of file ConversionProducer.h.

Referenced by ConversionProducer(), and preselectTrackPair().

double ConversionProducer::energyBC_
private

Definition at line 109 of file ConversionProducer.h.

Referenced by buildSuperAndBasicClusterGeoMap(), and ConversionProducer().

double ConversionProducer::energyTotalBC_
private

Definition at line 110 of file ConversionProducer.h.

Referenced by buildCollection(), checkTrackPair(), and ConversionProducer().

double ConversionProducer::halfWayEta_
private

Definition at line 102 of file ConversionProducer.h.

Referenced by ConversionProducer(), and getMatchedBC().

double ConversionProducer::halfWayPhi_
private

Definition at line 102 of file ConversionProducer.h.

double ConversionProducer::maxChi2Left_
private

Definition at line 115 of file ConversionProducer.h.

Referenced by ConversionProducer(), and trackQualityFilter().

double ConversionProducer::maxChi2Right_
private

Definition at line 115 of file ConversionProducer.h.

Referenced by ConversionProducer(), and trackQualityFilter().

unsigned int ConversionProducer::maxNumOfTrackInPU_
private

Definition at line 103 of file ConversionProducer.h.

Referenced by ConversionProducer(), and produce().

double ConversionProducer::maxTrackRho_
private

Definition at line 105 of file ConversionProducer.h.

Referenced by ConversionProducer(), and preselectTrackPair().

double ConversionProducer::maxTrackZ_
private

Definition at line 104 of file ConversionProducer.h.

Referenced by ConversionProducer(), and preselectTrackPair().

double ConversionProducer::minApproachHigh_
private

Definition at line 118 of file ConversionProducer.h.

Referenced by ConversionProducer(), and preselectTrackPair().

double ConversionProducer::minApproachLow_
private

Definition at line 118 of file ConversionProducer.h.

Referenced by ConversionProducer().

double ConversionProducer::minHitsLeft_
private

Definition at line 116 of file ConversionProducer.h.

Referenced by ConversionProducer(), and trackQualityFilter().

double ConversionProducer::minHitsRight_
private

Definition at line 116 of file ConversionProducer.h.

Referenced by ConversionProducer(), and trackQualityFilter().

double ConversionProducer::minSCEt_
private

Definition at line 106 of file ConversionProducer.h.

Referenced by buildSuperAndBasicClusterGeoMap(), and ConversionProducer().

double ConversionProducer::r_cut
private

Definition at line 121 of file ConversionProducer.h.

Referenced by ConversionProducer(), and preselectTrackPair().

bool ConversionProducer::rightBC_
private

Definition at line 125 of file ConversionProducer.h.

Referenced by checkTrackPair(), and ConversionProducer().

edm::EDGetTokenT<edm::View<reco::CaloCluster> > ConversionProducer::scBarrelProducer_
private

Definition at line 84 of file ConversionProducer.h.

Referenced by buildSuperAndBasicClusterGeoMap(), and ConversionProducer().

edm::EDGetTokenT<edm::View<reco::CaloCluster> > ConversionProducer::scEndcapProducer_
private

Definition at line 85 of file ConversionProducer.h.

Referenced by buildSuperAndBasicClusterGeoMap(), and ConversionProducer().

edm::EDGetTokenT<edm::View<reco::ConversionTrack> > ConversionProducer::src_
private

Definition at line 82 of file ConversionProducer.h.

Referenced by ConversionProducer(), and produce().

const TransientTrackBuilder* ConversionProducer::thettbuilder_
private

Definition at line 98 of file ConversionProducer.h.

Referenced by buildCollection(), ConversionProducer(), and produce().

ConversionVertexFinder* ConversionProducer::theVertexFinder_
private

Definition at line 96 of file ConversionProducer.h.

Referenced by checkVertex(), ConversionProducer(), and ~ConversionProducer().

bool ConversionProducer::usePvtx_
private

Definition at line 94 of file ConversionProducer.h.

Referenced by ConversionProducer(), and produce().

edm::EDGetTokenT<reco::VertexCollection> ConversionProducer::vertexProducer_
private

Definition at line 95 of file ConversionProducer.h.

Referenced by ConversionProducer(), and produce().

double ConversionProducer::vtxChi2_
private

Definition at line 122 of file ConversionProducer.h.

Referenced by buildCollection(), and ConversionProducer().