CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 &sumes)
 
void importToBlock (const edm::Event &, ElementList &) const override
 
- Public Member Functions inherited from BlockElementImporterBase
 BlockElementImporterBase (const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
 
 BlockElementImporterBase (const BlockElementImporterBase &)=delete
 
const std::string & name () const
 
BlockElementImporterBaseoperator= (const BlockElementImporterBase &)=delete
 
virtual void updateEventSetup (const edm::EventSetup &)
 

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
 

Detailed Description

Definition at line 14 of file GSFTrackImporter.cc.

Constructor & Destructor Documentation

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

Definition at line 16 of file GSFTrackImporter.cc.

17  :
18  BlockElementImporterBase(conf,sumes),
20  _isSecondary(conf.getParameter<bool>("gsfsAreSecondary")),
21  _superClustersArePF(conf.getParameter<bool>("superClustersArePF")){}
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
const bool _isSecondary
std::vector< GsfPFRecTrack > GsfPFRecTrackCollection
collection of GsfPFRecTrack objects
edm::EDGetTokenT< reco::GsfPFRecTrackCollection > _src
BlockElementImporterBase(const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
const bool _superClustersArePF

Member Function Documentation

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

Implements BlockElementImporterBase.

Definition at line 36 of file GSFTrackImporter.cc.

References _isSecondary, _src, _superClustersArePF, a, HLT_25ns14e33_v1_cff::distance, edm::Event::getByToken(), edm::Ref< C, T, F >::isAvailable(), edm::Ref< C, T, F >::isNonnull(), reco::PFBlockElement::SC, reco::PFBlockElementSuperCluster::setFromGsfElectron(), reco::PFBlockElementSuperCluster::setFromPFSuperCluster(), reco::PFBlockElementGsfTrack::setTrackType(), reco::PFBlockElement::T_FROM_GAMMACONV, and groupFilesInBlocks::temp.

37  {
40  e.getByToken(_src,gsftracks);
41  elems.reserve(elems.size() + gsftracks->size());
42  // setup our elements so that all the SCs are grouped together
43  auto SCs_end = std::partition(elems.begin(),elems.end(),
44  [](const ElementType& a){
45  return a->type() == reco::PFBlockElement::SC;
46  });
47  size_t SCs_end_position = std::distance(elems.begin(),SCs_end);
48  // insert gsf tracks and SCs, binding pre-existing SCs to ECAL-Driven GSF
49  auto bgsf = gsftracks->cbegin();
50  auto egsf = gsftracks->cend();
51  for( auto gsftrack = bgsf; gsftrack != egsf; ++gsftrack ) {
52  reco::GsfPFRecTrackRef gsfref(gsftracks,std::distance(bgsf,gsftrack));
53  const reco::GsfTrackRef& basegsfref = gsftrack->gsfTrackRef();
54  const auto& gsfextraref = basegsfref->extra();
55  // import associated super clusters
56  if( gsfextraref.isAvailable() && gsfextraref->seedRef().isAvailable()) {
57  reco::ElectronSeedRef seedref =
58  gsfextraref->seedRef().castTo<reco::ElectronSeedRef>();
59  if( seedref.isAvailable() && seedref->isEcalDriven() ) {
60  reco::SuperClusterRef scref =
61  seedref->caloCluster().castTo<reco::SuperClusterRef>();
62  if( scref.isNonnull() ) {
63  PFBlockElementSCEqual myEqual(scref);
64  auto sc_elem = std::find_if(elems.begin(),SCs_end,myEqual);
65  if( sc_elem != SCs_end ) {
67  static_cast<reco::PFBlockElementSuperCluster*>(sc_elem->get());
68  scbe->setFromGsfElectron(true);
69  } else {
72  scbe->setFromGsfElectron(true);
74  SCs_end = elems.insert(SCs_end,ElementType(scbe));
75  ++SCs_end; // point to element *after* the new one
76  }
77  }
78  }
79  }// gsf extra ref?
80  // cache the SC_end offset
81  SCs_end_position = std::distance(elems.begin(),SCs_end);
82  // get track momentum information
83  const std::vector<reco::PFTrajectoryPoint>& PfGsfPoint
84  = gsftrack->trajectoryPoints();
85  unsigned int c_gsf=0;
86  bool PassTracker = false;
87  bool GetPout = false;
88  unsigned int IndexPout = 0;
89  for(auto itPfGsfPoint = PfGsfPoint.begin();
90  itPfGsfPoint!= PfGsfPoint.end();++itPfGsfPoint) {
91  if (itPfGsfPoint->isValid()){
92  int layGsfP = itPfGsfPoint->layer();
93  if (layGsfP == -1) PassTracker = true;
94  if (PassTracker && layGsfP > 0 && GetPout == false) {
95  IndexPout = c_gsf-1;
96  GetPout = true;
97  }
98  //const math::XYZTLorentzVector GsfMoment = itPfGsfPoint->momentum();
99  ++c_gsf;
100  }
101  }
102  const math::XYZTLorentzVector& pin = PfGsfPoint[0].momentum();
103  const math::XYZTLorentzVector& pout = PfGsfPoint[IndexPout].momentum();
105  new reco::PFBlockElementGsfTrack(gsfref,pin,pout);
106  if( _isSecondary ) {
108  }
109  elems.emplace_back(temp);
110  // import brems from this primary gsf
111  for( const auto& brem : gsfref->PFRecBrem() ) {
112  const unsigned TrajP = brem.indTrajPoint();
113  if( TrajP != 99 ) {
114  const double DP = brem.DeltaP();
115  const double sDP = brem.SigmaDeltaP();
116  elems.emplace_back(new reco::PFBlockElementBrem(gsfref,DP,sDP,TrajP));
117  }
118  }
119  // protect against reallocations, create a fresh iterator
120  SCs_end = elems.begin() + SCs_end_position;
121  }// loop on gsf tracks
122  elems.shrink_to_fit();
123 }
bool isAvailable() const
Definition: Ref.h:576
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:252
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
const bool _isSecondary
void setFromGsfElectron(bool val)
set provenance
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
virtual void setTrackType(TrackType trType, bool value)
the trackType
Container::value_type value_type
double a
Definition: hdecay.h:121
edm::EDGetTokenT< reco::GsfPFRecTrackCollection > _src
const bool _superClustersArePF

Member Data Documentation

const bool GSFTrackImporter::_isSecondary
private

Definition at line 28 of file GSFTrackImporter.cc.

Referenced by importToBlock().

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

Definition at line 27 of file GSFTrackImporter.cc.

Referenced by importToBlock().

const bool GSFTrackImporter::_superClustersArePF
private

Definition at line 28 of file GSFTrackImporter.cc.

Referenced by importToBlock().