CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
GSFTrackImporter Class Reference
Inheritance diagram for GSFTrackImporter:
BlockElementImporterBase

Public Member Functions

 GSFTrackImporter (const edm::ParameterSet &conf, edm::ConsumesCollector &cc)
 
void importToBlock (const edm::Event &, ElementList &) const override
 
- Public Member Functions inherited from BlockElementImporterBase
 BlockElementImporterBase (const edm::ParameterSet &conf, edm::ConsumesCollector &cc)
 
 BlockElementImporterBase (const BlockElementImporterBase &)=delete
 
const std::string & name () const
 
BlockElementImporterBaseoperator= (const BlockElementImporterBase &)=delete
 
virtual void updateEventSetup (const edm::EventSetup &)
 
virtual ~BlockElementImporterBase ()=default
 

Private Attributes

const bool _isSecondary
 
edm::EDGetTokenT< reco::GsfPFRecTrackCollection_src
 
const bool _superClustersArePF
 

Additional Inherited Members

- Public Types inherited from BlockElementImporterBase
typedef std::vector< std::unique_ptr< reco::PFBlockElement > > ElementList
 
enum  VetoMode { pfRecTrackCollection = 0, ticlSeedingRegion = 1, pfCandidateCollection = 2 }
 

Detailed Description

Definition at line 11 of file GSFTrackImporter.cc.

Constructor & Destructor Documentation

◆ GSFTrackImporter()

GSFTrackImporter::GSFTrackImporter ( const edm::ParameterSet conf,
edm::ConsumesCollector cc 
)
inline

Definition at line 13 of file GSFTrackImporter.cc.

16  _isSecondary(conf.getParameter<bool>("gsfsAreSecondary")),
17  _superClustersArePF(conf.getParameter<bool>("superClustersArePF")) {}
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
const bool _isSecondary
BlockElementImporterBase(const edm::ParameterSet &conf, edm::ConsumesCollector &cc)
std::vector< GsfPFRecTrack > GsfPFRecTrackCollection
collection of GsfPFRecTrack objects
edm::EDGetTokenT< reco::GsfPFRecTrackCollection > _src
const bool _superClustersArePF

Member Function Documentation

◆ importToBlock()

void GSFTrackImporter::importToBlock ( const edm::Event e,
BlockElementImporterBase::ElementList elems 
) const
overridevirtual

Implements BlockElementImporterBase.

Definition at line 28 of file GSFTrackImporter.cc.

References _isSecondary, _src, _superClustersArePF, a, HLT_2024v14_cff::distance, MillePedeFileConverter_cfg::e, DetId::Forward, DetId::HGCalEE, DetId::HGCalHSc, DetId::HGCalHSi, edm::Ref< C, T, F >::isAvailable(), edm::Ref< C, T, F >::isNonnull(), l1ctLayer1_patternWriters_cff::partition, reco::PFBlockElement::SC, reco::PFBlockElementSuperCluster::setFromGsfElectron(), reco::PFBlockElementSuperCluster::setFromPFSuperCluster(), reco::PFBlockElement::T_FROM_GAMMACONV, and groupFilesInBlocks::temp.

28  {
29  auto gsftracks = e.getHandle(_src);
30  elems.reserve(elems.size() + gsftracks->size());
31  // setup our elements so that all the SCs are grouped together
32  auto SCs_end =
33  std::partition(elems.begin(), elems.end(), [](const auto& a) { return a->type() == reco::PFBlockElement::SC; });
34  // insert gsf tracks and SCs, binding pre-existing SCs to ECAL-Driven GSF
35  auto bgsf = gsftracks->cbegin();
36  auto egsf = gsftracks->cend();
37  for (auto gsftrack = bgsf; gsftrack != egsf; ++gsftrack) {
38  reco::GsfPFRecTrackRef gsfref(gsftracks, std::distance(bgsf, gsftrack));
39  const reco::GsfTrackRef& basegsfref = gsftrack->gsfTrackRef();
40  const auto& gsfextraref = basegsfref->extra();
41  // import associated super clusters
42  if (gsfextraref.isAvailable() && gsfextraref->seedRef().isAvailable()) {
43  reco::ElectronSeedRef seedref = gsfextraref->seedRef().castTo<reco::ElectronSeedRef>();
44  if (seedref.isAvailable() && seedref->isEcalDriven()) {
45  reco::SuperClusterRef scref = seedref->caloCluster().castTo<reco::SuperClusterRef>();
46  if (scref.isNonnull()) {
47  // explicitly veto HGCal super clusters
48  if (scref->seed()->seed().det() == DetId::HGCalEE || scref->seed()->seed().det() == DetId::HGCalHSi ||
49  scref->seed()->seed().det() == DetId::HGCalHSc || scref->seed()->seed().det() == DetId::Forward)
50  continue;
51  PFBlockElementSCEqual myEqual(scref);
52  auto sc_elem = std::find_if(elems.begin(), SCs_end, myEqual);
53  if (sc_elem != SCs_end) {
54  reco::PFBlockElementSuperCluster* scbe = static_cast<reco::PFBlockElementSuperCluster*>(sc_elem->get());
55  scbe->setFromGsfElectron(true);
56  } else {
58  scbe->setFromGsfElectron(true);
60  SCs_end = elems.emplace(SCs_end, scbe);
61  ++SCs_end; // point to element *after* the new one
62  }
63  }
64  }
65  } // gsf extra ref?
66  // cache the SC_end offset
67  size_t SCs_end_position = std::distance(elems.begin(), SCs_end);
68  // get track momentum information
69  const std::vector<reco::PFTrajectoryPoint>& PfGsfPoint = gsftrack->trajectoryPoints();
70  unsigned int c_gsf = 0;
71  bool PassTracker = false;
72  unsigned int IndexPout = 0;
73  for (const auto& pfGsfPoint : PfGsfPoint) {
74  if (pfGsfPoint.isValid()) {
75  int layGsfP = pfGsfPoint.layer();
76  if (layGsfP == -1)
77  PassTracker = true;
78  else if (PassTracker && layGsfP > 0) {
79  IndexPout = c_gsf - 1;
80  break;
81  }
82  ++c_gsf;
83  }
84  }
85  const math::XYZTLorentzVector& pin = PfGsfPoint[0].momentum();
86  const math::XYZTLorentzVector& pout = PfGsfPoint[IndexPout].momentum();
88  if (_isSecondary) {
89  temp->setTrackType(reco::PFBlockElement::T_FROM_GAMMACONV, true);
90  }
91  elems.emplace_back(temp);
92  // import brems from this primary gsf
93  for (const auto& brem : gsfref->PFRecBrem()) {
94  const unsigned TrajP = brem.indTrajPoint();
95  if (TrajP != 99) {
96  const double DP = brem.DeltaP();
97  const double sDP = brem.SigmaDeltaP();
98  elems.emplace_back(new reco::PFBlockElementBrem(gsfref, DP, sDP, TrajP));
99  }
100  }
101  // protect against reallocations, create a fresh iterator
102  SCs_end = elems.begin() + SCs_end_position;
103  } // loop on gsf tracks
104  elems.shrink_to_fit();
105 }
const bool _isSecondary
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:232
void setFromGsfElectron(bool val)
set provenance
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
bool isAvailable() const
Definition: Ref.h:535
double a
Definition: hdecay.h:121
edm::EDGetTokenT< reco::GsfPFRecTrackCollection > _src
const bool _superClustersArePF

Member Data Documentation

◆ _isSecondary

const bool GSFTrackImporter::_isSecondary
private

Definition at line 23 of file GSFTrackImporter.cc.

Referenced by importToBlock().

◆ _src

edm::EDGetTokenT<reco::GsfPFRecTrackCollection> GSFTrackImporter::_src
private

Definition at line 22 of file GSFTrackImporter.cc.

Referenced by importToBlock().

◆ _superClustersArePF

const bool GSFTrackImporter::_superClustersArePF
private

Definition at line 23 of file GSFTrackImporter.cc.

Referenced by importToBlock().