CMS 3D CMS Logo

GSFTrackImporter.cc
Go to the documentation of this file.
13 
15 public:
17  edm::ConsumesCollector& sumes) :
18  BlockElementImporterBase(conf,sumes),
19  _src(sumes.consumes<reco::GsfPFRecTrackCollection>(conf.getParameter<edm::InputTag>("source"))),
20  _isSecondary(conf.getParameter<bool>("gsfsAreSecondary")),
21  _superClustersArePF(conf.getParameter<bool>("superClustersArePF")){}
22 
23  void importToBlock( const edm::Event& ,
24  ElementList& ) const override;
25 
26 private:
29 };
30 
33  "GSFTrackImporter");
34 
38  auto gsftracks = e.getHandle(_src);
39  elems.reserve(elems.size() + gsftracks->size());
40  // setup our elements so that all the SCs are grouped together
41  auto SCs_end = std::partition(elems.begin(),elems.end(),
42  [](const auto& a){
43  return a->type() == reco::PFBlockElement::SC;
44  });
45  // insert gsf tracks and SCs, binding pre-existing SCs to ECAL-Driven GSF
46  auto bgsf = gsftracks->cbegin();
47  auto egsf = gsftracks->cend();
48  for( auto gsftrack = bgsf; gsftrack != egsf; ++gsftrack ) {
49  reco::GsfPFRecTrackRef gsfref(gsftracks,std::distance(bgsf,gsftrack));
50  const reco::GsfTrackRef& basegsfref = gsftrack->gsfTrackRef();
51  const auto& gsfextraref = basegsfref->extra();
52  // import associated super clusters
53  if( gsfextraref.isAvailable() && gsfextraref->seedRef().isAvailable()) {
54  reco::ElectronSeedRef seedref =
55  gsfextraref->seedRef().castTo<reco::ElectronSeedRef>();
56  if( seedref.isAvailable() && seedref->isEcalDriven() ) {
57  reco::SuperClusterRef scref =
58  seedref->caloCluster().castTo<reco::SuperClusterRef>();
59  if( scref.isNonnull() ) {
60  // explicitly veto HGCal super clusters
61  if( scref->seed()->seed().det() == DetId::Forward ) continue;
62  PFBlockElementSCEqual myEqual(scref);
63  auto sc_elem = std::find_if(elems.begin(),SCs_end,myEqual);
64  if( sc_elem != SCs_end ) {
66  static_cast<reco::PFBlockElementSuperCluster*>(sc_elem->get());
67  scbe->setFromGsfElectron(true);
68  } else {
71  scbe->setFromGsfElectron(true);
73  SCs_end = elems.emplace(SCs_end,scbe);
74  ++SCs_end; // point to element *after* the new one
75  }
76  }
77  }
78  }// gsf extra ref?
79  // cache the SC_end offset
80  size_t SCs_end_position = std::distance(elems.begin(),SCs_end);
81  // get track momentum information
82  const std::vector<reco::PFTrajectoryPoint>& PfGsfPoint = gsftrack->trajectoryPoints();
83  unsigned int c_gsf=0;
84  bool PassTracker = false;
85  unsigned int IndexPout = 0;
86  for(const auto& pfGsfPoint : PfGsfPoint) {
87  if (pfGsfPoint.isValid()){
88  int layGsfP = pfGsfPoint.layer();
89  if (layGsfP == -1) PassTracker = true;
90  else if (PassTracker && layGsfP > 0) {
91  IndexPout = c_gsf-1;
92  break;
93  }
94  ++c_gsf;
95  }
96  }
97  const math::XYZTLorentzVector& pin = PfGsfPoint[0].momentum();
98  const math::XYZTLorentzVector& pout = PfGsfPoint[IndexPout].momentum();
100  new reco::PFBlockElementGsfTrack(gsfref,pin,pout);
101  if( _isSecondary ) {
103  }
104  elems.emplace_back(temp);
105  // import brems from this primary gsf
106  for( const auto& brem : gsfref->PFRecBrem() ) {
107  const unsigned TrajP = brem.indTrajPoint();
108  if( TrajP != 99 ) {
109  const double DP = brem.DeltaP();
110  const double sDP = brem.SigmaDeltaP();
111  elems.emplace_back(new reco::PFBlockElementBrem(gsfref,DP,sDP,TrajP));
112  }
113  }
114  // protect against reallocations, create a fresh iterator
115  SCs_end = elems.begin() + SCs_end_position;
116  }// loop on gsf tracks
117  elems.shrink_to_fit();
118 }
bool isAvailable() const
Definition: Ref.h:575
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:251
const bool _isSecondary
void setFromGsfElectron(bool val)
set provenance
GSFTrackImporter(const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
Handle< PROD > getHandle(EDGetTokenT< PROD > token) const
Definition: Event.h:539
std::vector< GsfPFRecTrack > GsfPFRecTrackCollection
collection of GsfPFRecTrack objects
void importToBlock(const edm::Event &, ElementList &) const override
void setTrackType(TrackType trType, bool value) override
the trackType
fixed size matrix
HLT enums.
double a
Definition: hdecay.h:121
#define DEFINE_EDM_PLUGIN(factory, type, name)
std::vector< std::unique_ptr< reco::PFBlockElement > > ElementList
edm::EDGetTokenT< reco::GsfPFRecTrackCollection > _src
const bool _superClustersArePF