CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
LowPtGsfElectronSCProducer Class Reference

#include <LowPtGsfElectronSCProducer.h>

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

Public Member Functions

 LowPtGsfElectronSCProducer (const edm::ParameterSet &)
 
void produce (edm::Event &, const edm::EventSetup &) override
 
 ~LowPtGsfElectronSCProducer () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
bool hasAbilityToProduceInLumis () const final
 
bool hasAbilityToProduceInRuns () const final
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &)
 

Private Member Functions

reco::PFClusterRef closestCluster (const reco::PFTrajectoryPoint &point, const edm::Handle< reco::PFClusterCollection > &clusters, std::vector< int > &matched)
 

Private Attributes

const double dr2_
 
const edm::EDGetTokenT< reco::PFClusterCollectionecalClusters_
 
const edm::EDGetTokenT< reco::GsfPFRecTrackCollectiongsfPfRecTracks_
 

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

Definition at line 13 of file LowPtGsfElectronSCProducer.h.

Constructor & Destructor Documentation

LowPtGsfElectronSCProducer::LowPtGsfElectronSCProducer ( const edm::ParameterSet cfg)
explicit

Definition at line 20 of file LowPtGsfElectronSCProducer.cc.

References looper::cfg, dr2_, ecalClusters_, and edm::ParameterSet::getParameter().

20  :
21  gsfPfRecTracks_{consumes<reco::GsfPFRecTrackCollection>( cfg.getParameter<edm::InputTag>("gsfPfRecTracks") )},
22  ecalClusters_{consumes<reco::PFClusterCollection>( cfg.getParameter<edm::InputTag>("ecalClusters") )},
23  dr2_{cfg.getParameter<double>("MaxDeltaR2")}
24 {
25  produces<reco::CaloClusterCollection>();
26  produces<reco::SuperClusterCollection>();
27  produces< edm::ValueMap<reco::SuperClusterRef> >();
28 }
T getParameter(std::string const &) const
const edm::EDGetTokenT< reco::PFClusterCollection > ecalClusters_
const edm::EDGetTokenT< reco::GsfPFRecTrackCollection > gsfPfRecTracks_
LowPtGsfElectronSCProducer::~LowPtGsfElectronSCProducer ( )
override

Definition at line 32 of file LowPtGsfElectronSCProducer.cc.

33 {}

Member Function Documentation

reco::PFClusterRef LowPtGsfElectronSCProducer::closestCluster ( const reco::PFTrajectoryPoint point,
const edm::Handle< reco::PFClusterCollection > &  clusters,
std::vector< int > &  matched 
)
private

Definition at line 200 of file LowPtGsfElectronSCProducer.cc.

References reco::deltaR2(), dr2_, MillePedeFileConverter_cfg::e, spr::find(), cuy::ii, edm::Ref< C, T, F >::index(), reco::PFTrajectoryPoint::isValid(), and reco::PFTrajectoryPoint::positionREP().

202  {
203  reco::PFClusterRef closest;
204  if ( point.isValid() ) {
205  float dr2min = dr2_;
206  for ( size_t ii = 0; ii < clusters->size(); ++ii ) {
207  if ( std::find( matched.begin(), matched.end(), ii ) == matched.end() ) {
208  float dr2 = reco::deltaR2( clusters->at(ii), point.positionREP() );
209  if ( dr2 < dr2min ) {
210  closest = reco::PFClusterRef( clusters, ii );
211  dr2min = dr2;
212  }
213  }
214  }
215  if ( dr2min < (dr2_ - 1.e-6) ) { matched.push_back( closest.index() ); }
216  }
217  return closest;
218 }
const REPPoint & positionREP() const
trajectory position in (rho, eta, phi) base
key_type index() const
Definition: Ref.h:268
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
edm::Ref< PFClusterCollection > PFClusterRef
persistent reference to PFCluster objects
Definition: PFClusterFwd.h:15
ii
Definition: cuy.py:589
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
bool isValid() const
is this point valid ?
void LowPtGsfElectronSCProducer::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 222 of file LowPtGsfElectronSCProducer.cc.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), and DEFINE_FWK_MODULE.

223 {
225  desc.add<edm::InputTag>("gsfPfRecTracks",edm::InputTag("lowPtGsfElePfGsfTracks"));
226  desc.add<edm::InputTag>("ecalClusters",edm::InputTag("particleFlowClusterECAL"));
227  desc.add<edm::InputTag>("hcalClusters",edm::InputTag("particleFlowClusterHCAL"));
228  desc.add<double>("MaxDeltaR2",0.5);
229  descriptions.add("defaultLowPtGsfElectronSuperClusters",desc);
230 }
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void LowPtGsfElectronSCProducer::produce ( edm::Event event,
const edm::EventSetup setup 
)
override

Definition at line 37 of file LowPtGsfElectronSCProducer.cc.

References fastPrimaryVertexProducer_cfi::clusters, reco::deltaR2(), dr2_, RecoEcal_cff::ecalClusters, ecalClusters_, edm::helper::Filler< Map >::fill(), objects.autophobj::filler, gsfElectronCores_cfi::gsfPfRecTracks, gsfPfRecTracks_, edm::RefProd< C >::id(), edm::helper::Filler< Map >::insert(), plotBeamSpotDB::ipoint, edm::Ptr< T >::isNonnull(), edm::Ptr< T >::isNull(), edm::HandleBase::isValid(), eostools::move(), PFClusterWidthAlgo::pflowEtaWidth(), PFClusterWidthAlgo::pflowPhiWidth(), point, hiPixelPairStep_cff::points, edm::PtrVector< T >::push_back(), reco::SuperCluster::rawEnergy(), SimDataFormats::CaloAnalysis::sc, SurveyInfoScenario_cff::seed, reco::SuperCluster::setClusters(), reco::CaloCluster::setCorrectedEnergy(), reco::SuperCluster::setEtaWidth(), reco::SuperCluster::setPhiWidth(), reco::SuperCluster::setSeed(), findQualityFiles::size, X, DOFs::Y, and DOFs::Z.

38 {
39 
40  // Input GsfPFRecTracks collection
42  event.getByToken(gsfPfRecTracks_,gsfPfRecTracks);
43  if ( !gsfPfRecTracks.isValid() ) { edm::LogError("Problem with gsfPfRecTracks handle"); }
44 
45  // Input EcalClusters collection
47  event.getByToken(ecalClusters_,ecalClusters);
48  if ( !ecalClusters.isValid() ) { edm::LogError("Problem with ecalClusters handle"); }
49 
50  // Output SuperClusters collection and getRefBeforePut
51  auto superClusters = std::make_unique<reco::SuperClusterCollection>( reco::SuperClusterCollection() );
52  superClusters->reserve(gsfPfRecTracks->size());
53  const reco::SuperClusterRefProd superClustersRefProd = event.getRefBeforePut<reco::SuperClusterCollection>();
54 
55  // Output ValueMap container of GsfPFRecTrackRef index to SuperClusterRef
56  std::vector<reco::SuperClusterRef> superClustersValueMap;
57 
58  // Output CaloClusters collection
59  auto caloClusters = std::make_unique<reco::CaloClusterCollection>( reco::CaloClusterCollection() );
60  caloClusters->reserve(ecalClusters->size());
61 
62  // Index[GSF track][trajectory point] for "best" CaloCluster
63  std::vector< std::vector<int> > cluster_idx;
64  cluster_idx.reserve( gsfPfRecTracks->size());
65 
66  // Index[GSF track][trajectory point] for "best" PFCluster
67  std::vector< std::vector<int> > pfcluster_idx;
68  pfcluster_idx.reserve( gsfPfRecTracks->size());
69 
70  // dr2min[GSF track][trajectory point] for "best" CaloCluster
71  std::vector< std::vector<float> > cluster_dr2min;
72  cluster_dr2min.reserve( gsfPfRecTracks->size());
73 
74  // Construct list of trajectory points from the GSF track and electron brems
75  std::vector< std::vector<const reco::PFTrajectoryPoint*> > points;
76  points.reserve( gsfPfRecTracks->size());
77  for ( auto const& trk: *gsfPfRecTracks) {
78  // Extrapolated track
79  std::vector<const reco::PFTrajectoryPoint*> traj;
80  traj.reserve(trk.PFRecBrem().size()+1);
81  traj.push_back( &trk.extrapolatedPoint(reco::PFTrajectoryPoint::LayerType::ECALShowerMax) );
82  // Extrapolated brem trajectories
83  for ( auto const& brem : trk.PFRecBrem() ) {
84  traj.push_back( &brem.extrapolatedPoint(reco::PFTrajectoryPoint::LayerType::ECALShowerMax) );
85  }
86  auto size = traj.size();
87  points.push_back(std::move(traj));
88  // Size containers
89  cluster_idx.emplace_back(size,-1);
90  pfcluster_idx.emplace_back(size,-1);
91  cluster_dr2min.emplace_back(size,1.e6);
92  }
93 
94  // For each cluster, find closest trajectory point ("global" arbitration at event level)
95  for ( size_t iclu = 0; iclu < ecalClusters->size(); ++iclu ) { // Cluster loop
96  std::pair<int,int> point = std::make_pair(-1,-1);
97  float dr2min = 1.e6;
98  for ( size_t ipoint = 0; ipoint < points.size(); ++ipoint ) { // GSF track loop
99  for ( size_t jpoint = 0; jpoint < points[ipoint].size(); ++jpoint ) { // Traj point loop
100  if ( points[ipoint][jpoint]->isValid() ) {
101  float dr2 = reco::deltaR2( ecalClusters->at(iclu), points[ipoint][jpoint]->positionREP() );
102  if ( dr2 < dr2min ) {
103  // Store nearest point to this cluster
104  dr2min = dr2;
105  point = std::make_pair(ipoint,jpoint);
106  }
107  }
108  }
109  }
110  if ( point.first >= 0 && point.second >= 0 && // if this cluster is matched to a point ...
111  dr2min < cluster_dr2min[point.first][point.second] ) { // ... and cluster is closest to the same point
112  // Copy CaloCluster to new collection
113  caloClusters->push_back(ecalClusters->at(iclu));
114  // Store cluster index for creation of SC later
115  cluster_idx[point.first][point.second] = caloClusters->size()-1;
116  pfcluster_idx[point.first][point.second] = iclu;
117  cluster_dr2min[point.first][point.second] = dr2min;
118  }
119  }
120 
121  // Put CaloClusters in event and get orphan handle
122  const edm::OrphanHandle<reco::CaloClusterCollection>& caloClustersH = event.put(std::move(caloClusters));
123 
124  // Loop through GSF tracks
125  for ( size_t itrk = 0; itrk < gsfPfRecTracks->size(); ++itrk ) {
126 
127  // Used to create SC
128  float energy = 0.;
129  float X = 0., Y = 0., Z = 0.;
132  std::vector<const reco::PFCluster*> barePtrs;
133 
134  // Find closest match in dr2 from points associated to given track
135  int index = -1;
136  float dr2 = 1.e6;
137  for ( size_t ipoint = 0; ipoint < cluster_idx[itrk].size(); ++ipoint ) {
138  if ( cluster_idx[itrk][ipoint] < 0 ) { continue; }
139  if ( cluster_dr2min[itrk][ipoint] < dr2 ) {
140  dr2 = cluster_dr2min[itrk][ipoint];
141  index = cluster_idx[itrk][ipoint];
142  }
143  }
144 
145  // For each track, loop through points and use associated cluster
146  for ( size_t ipoint = 0; ipoint < cluster_idx[itrk].size(); ++ipoint ) {
147  if ( cluster_idx[itrk][ipoint] < 0 ) { continue; }
148  reco::CaloClusterPtr clu(caloClustersH,cluster_idx[itrk][ipoint]);
149  if ( clu.isNull() ) { continue; }
150  if ( !( cluster_dr2min[itrk][ipoint] < dr2_ || // Require cluster to be closer than dr2_ ...
151  index == cluster_idx[itrk][ipoint] ) ) { continue; } // ... unless it is the closest one ...
152  if ( seed.isNull() ) { seed = clu; }
153  clusters.push_back(clu);
154  energy += clu->correctedEnergy();
155  X += clu->position().X() * clu->correctedEnergy();
156  Y += clu->position().Y() * clu->correctedEnergy();
157  Z += clu->position().Z() * clu->correctedEnergy();
158  auto index = pfcluster_idx[itrk][ipoint];
159  if(index < static_cast<decltype(index)>(ecalClusters->size())) {
160  barePtrs.push_back(&(*ecalClusters)[index]);
161  }
162  }
163  X /= energy;
164  Y /= energy;
165  Z /= energy;
166 
167  // Create SC
168  if ( seed.isNonnull() ) {
170  sc.setCorrectedEnergy(energy);
171  sc.setSeed(seed);
172  sc.setClusters(clusters);
173  PFClusterWidthAlgo pfwidth(barePtrs);
174  sc.setEtaWidth(pfwidth.pflowEtaWidth());
175  sc.setPhiWidth(pfwidth.pflowPhiWidth());
176  sc.rawEnergy(); // Cache the value of raw energy
177  superClusters->push_back(sc);
178 
179  // Populate ValueMap container
180  superClustersValueMap.push_back( reco::SuperClusterRef( superClustersRefProd, superClusters->size()-1 ) );
181  } else {
182  superClustersValueMap.push_back( reco::SuperClusterRef( superClustersRefProd.id() ) );
183  }
184 
185  } // GSF tracks
186 
187  // Put SuperClusters in event
188  event.put(std::move(superClusters));
189 
190  auto ptr = std::make_unique< edm::ValueMap<reco::SuperClusterRef> >( edm::ValueMap<reco::SuperClusterRef>() );
192  filler.insert(gsfPfRecTracks, superClustersValueMap.begin(), superClustersValueMap.end());
193  filler.fill();
194  event.put(std::move(ptr));
195 
196 }
size
Write out results.
const edm::EDGetTokenT< reco::PFClusterCollection > ecalClusters_
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:140
#define X(str)
Definition: MuonsGrabber.cc:48
const edm::EDGetTokenT< reco::GsfPFRecTrackCollection > gsfPfRecTracks_
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
bool isNull() const
Checks for null.
Definition: Ptr.h:164
std::vector< CaloCluster > CaloClusterCollection
collection of CaloCluster objects
bool isValid() const
Definition: HandleBase.h:74
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:168
ProductID id() const
Accessor for product ID.
Definition: RefProd.h:137
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
def move(src, dest)
Definition: eostools.py:510
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5

Member Data Documentation

const double LowPtGsfElectronSCProducer::dr2_
private
const edm::EDGetTokenT<reco::PFClusterCollection> LowPtGsfElectronSCProducer::ecalClusters_
private

Definition at line 32 of file LowPtGsfElectronSCProducer.h.

Referenced by LowPtGsfElectronSCProducer(), and produce().

const edm::EDGetTokenT<reco::GsfPFRecTrackCollection> LowPtGsfElectronSCProducer::gsfPfRecTracks_
private

Definition at line 31 of file LowPtGsfElectronSCProducer.h.

Referenced by produce().