CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
calib::CalibElectron Class Reference

#include <CalibElectron.h>

Public Member Functions

 CalibElectron ()
 
 CalibElectron (const reco::GsfElectron *ele, const EcalRecHitCollection *theHits, const EcalRecHitCollection *theEEHits)
 
std::vector< std::pair< int, float > > getCalibModulesWeights (TString calibtype)
 
const EcalRecHitCollectiongetEERecHits ()
 
const EcalRecHitCollectiongetRecHits ()
 
const reco::GsfElectrongetRecoElectron ()
 
 ~CalibElectron ()
 

Private Attributes

const EcalRecHitCollectiontheEEHits_
 
const reco::GsfElectrontheElectron_
 
const EcalRecHitCollectiontheHits_
 

Detailed Description

Definition at line 13 of file CalibElectron.h.

Constructor & Destructor Documentation

CalibElectron::CalibElectron ( )

Definition at line 14 of file CalibElectron.cc.

14 : theElectron_(nullptr), theHits_(nullptr), theEEHits_(nullptr) {}
const EcalRecHitCollection * theEEHits_
Definition: CalibElectron.h:32
const EcalRecHitCollection * theHits_
Definition: CalibElectron.h:31
const reco::GsfElectron * theElectron_
Definition: CalibElectron.h:29
calib::CalibElectron::CalibElectron ( const reco::GsfElectron ele,
const EcalRecHitCollection theHits,
const EcalRecHitCollection theEEHits 
)
inline

Definition at line 16 of file CalibElectron.h.

19  : theElectron_(ele), theHits_(theHits), theEEHits_(theEEHits){};
const EcalRecHitCollection * theEEHits_
Definition: CalibElectron.h:32
const EcalRecHitCollection * theHits_
Definition: CalibElectron.h:31
const reco::GsfElectron * theElectron_
Definition: CalibElectron.h:29
calib::CalibElectron::~CalibElectron ( )
inline

Definition at line 21 of file CalibElectron.h.

References getCalibModulesWeights().

21 {};

Member Function Documentation

std::vector< std::pair< int, float > > CalibElectron::getCalibModulesWeights ( TString  calibtype)

Definition at line 16 of file CalibElectron.cc.

References gather_cfg::cout, EcalBarrel, EcalEndcap, EcalRecHit::energy(), EgHLTOffHistBins_cfi::et, PVValHelper::eta, reco::LeafCandidate::eta(), JetChargeProducer_cfi::exp, edm::SortedCollection< T, SORT >::find(), EcalIndexingTools::getInstance(), EcalRingCalibrationTools::getModuleIndex(), EcalIndexingTools::getNumberOfChannels(), EcalIndexingTools::getProgressiveIndex(), EcalRingCalibrationTools::getRingIndex(), mps_fire::i, recoMuon::in, EcalRingCalibrationTools::N_MODULES_BARREL, EcalRingCalibrationTools::N_RING_TOTAL, funct::sin(), reco::GsfElectron::superCluster(), theEEHits_, theElectron_, theHits_, and theta().

Referenced by ~CalibElectron().

16  {
17  std::vector<std::pair<int, float> > theWeights;
18 
19  if (calibtype == "RING") {
21 
22  for (int i = 0; i < EcalRingCalibrationTools::N_RING_TOTAL; ++i)
23  w_ring[i] = 0.;
24 
25  std::vector<std::pair<DetId, float> > scDetIds = theElectron_->superCluster()->hitsAndFractions();
26 
27  for (std::vector<std::pair<DetId, float> >::const_iterator idIt = scDetIds.begin(); idIt != scDetIds.end();
28  ++idIt) {
29  const EcalRecHit* rh = nullptr;
30  if ((*idIt).first.subdetId() == EcalBarrel)
31  rh = &*(theHits_->find((*idIt).first));
32  else if ((*idIt).first.subdetId() == EcalEndcap)
33  rh = &*(theEEHits_->find((*idIt).first));
34  if (!rh)
35  std::cout << "CalibElectron::BIG ERROR::RecHit NOT FOUND" << std::endl;
36  w_ring[EcalRingCalibrationTools::getRingIndex((*idIt).first)] += rh->energy();
37  }
38 
39  for (int i = 0; i < EcalRingCalibrationTools::N_RING_TOTAL; ++i)
40  if (w_ring[i] != 0.)
41  theWeights.push_back(std::pair<int, float>(i, w_ring[i]));
42  // std::cout << " ring " << i << " - energy sum " << w_ring[i] << std::endl;
43  }
44 
45  else if (calibtype == "MODULE") {
47 
49  w_ring[i] = 0.;
50 
51  std::vector<std::pair<DetId, float> > scDetIds = theElectron_->superCluster()->hitsAndFractions();
52 
53  for (std::vector<std::pair<DetId, float> >::const_iterator idIt = scDetIds.begin(); idIt != scDetIds.end();
54  ++idIt) {
55  const EcalRecHit* rh = nullptr;
56  if ((*idIt).first.subdetId() == EcalBarrel)
57  rh = &*(theHits_->find((*idIt).first));
58  else if ((*idIt).first.subdetId() == EcalEndcap)
59  rh = &*(theEEHits_->find((*idIt).first));
60  if (!rh)
61  std::cout << "CalibElectron::BIG ERROR::RecHit NOT FOUND" << std::endl;
62  w_ring[EcalRingCalibrationTools::getModuleIndex((*idIt).first)] += rh->energy();
63  }
64 
66  if (w_ring[i] != 0.)
67  theWeights.push_back(std::pair<int, float>(i, w_ring[i]));
68  // std::cout << " ring " << i << " - energy sum " << w_ring[i] << std::endl;
69 
70  }
71 
72  else if (calibtype == "ABS_SCALE") {
73  std::cout << "ENTERING CalibElectron, ABS SCALE mode" << std::endl;
74 
75  float w_ring(0.);
76 
77  std::vector<std::pair<DetId, float> > scDetIds = theElectron_->superCluster()->hitsAndFractions();
78 
79  for (std::vector<std::pair<DetId, float> >::const_iterator idIt = scDetIds.begin(); idIt != scDetIds.end();
80  ++idIt) {
81  const EcalRecHit* rh = nullptr;
82  if ((*idIt).first.subdetId() == EcalBarrel)
83  rh = &*(theHits_->find((*idIt).first));
84  else if ((*idIt).first.subdetId() == EcalEndcap)
85  rh = &*(theEEHits_->find((*idIt).first));
86  if (!rh)
87  std::cout << "CalibElectron::BIG ERROR::RecHit NOT FOUND" << std::endl;
88 
89  w_ring += rh->energy();
90  }
91 
92  if (w_ring != 0.)
93  theWeights.push_back(std::pair<int, float>(0, w_ring));
94  std::cout << " ABS SCALE - energy sum " << w_ring << std::endl;
95 
96  }
97 
98  else if (calibtype == "ETA_ET_MODE") {
99  float w_ring[200];
100 
101  for (int i = 0; i < EcalIndexingTools::getInstance()->getNumberOfChannels(); ++i)
102  w_ring[i] = 0.;
103 
104  std::vector<std::pair<DetId, float> > scDetIds = theElectron_->superCluster()->hitsAndFractions();
105 
106  for (std::vector<std::pair<DetId, float> >::const_iterator idIt = scDetIds.begin(); idIt != scDetIds.end();
107  ++idIt) {
108  const EcalRecHit* rh = nullptr;
109  if ((*idIt).first.subdetId() == EcalBarrel)
110  rh = &*(theHits_->find((*idIt).first));
111  else if ((*idIt).first.subdetId() == EcalEndcap)
112  rh = &*(theEEHits_->find((*idIt).first));
113  if (!rh)
114  std::cout << "CalibElectron::BIG ERROR::RecHit NOT FOUND" << std::endl;
115 
116  float eta = fabs(theElectron_->eta());
117  float theta = 2. * atan(exp(-eta));
118  float et = theElectron_->superCluster()->energy() * sin(theta);
119 
121 
122  w_ring[in] += rh->energy();
123  //w_ring[in]+=theElectron_->superCluster()->energy();
124 
125  std::cout << "CalibElectron::filling channel " << in << " with value " << theElectron_->superCluster()->energy()
126  << std::endl;
127  }
128 
129  for (int i = 0; i < EcalIndexingTools::getInstance()->getNumberOfChannels(); ++i) {
130  if (w_ring[i] != 0.) {
131  theWeights.push_back(std::pair<int, float>(i, w_ring[i]));
132  std::cout << " ring " << i << " - energy sum " << w_ring[i] << std::endl;
133  }
134  }
135 
136  }
137 
138  else {
139  std::cout << "CalibType not yet implemented" << std::endl;
140  }
141 
142  return theWeights;
143 }
static short getModuleIndex(DetId aDetId)
double eta() const final
momentum pseudorapidity
const EcalRecHitCollection * theEEHits_
Definition: CalibElectron.h:32
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
static short getRingIndex(DetId aDetId)
Retrieve the phi-ring index corresponding to a DetId.
const EcalRecHitCollection * theHits_
Definition: CalibElectron.h:31
static EcalIndexingTools * getInstance()
float energy() const
Definition: EcalRecHit.h:68
int getProgressiveIndex(double, double)
iterator find(key_type k)
SuperClusterRef superCluster() const override
reference to a SuperCluster
Definition: GsfElectron.h:155
const reco::GsfElectron * theElectron_
Definition: CalibElectron.h:29
static constexpr short N_RING_TOTAL
static constexpr short N_MODULES_BARREL
const EcalRecHitCollection* calib::CalibElectron::getEERecHits ( )
inline

Definition at line 26 of file CalibElectron.h.

References theEEHits_.

26 { return theEEHits_; }
const EcalRecHitCollection * theEEHits_
Definition: CalibElectron.h:32
const EcalRecHitCollection* calib::CalibElectron::getRecHits ( )
inline

Definition at line 25 of file CalibElectron.h.

References theHits_.

25 { return theHits_; }
const EcalRecHitCollection * theHits_
Definition: CalibElectron.h:31
const reco::GsfElectron* calib::CalibElectron::getRecoElectron ( )
inline

Definition at line 24 of file CalibElectron.h.

References theElectron_.

Referenced by ZIterativeAlgorithmWithFit::addEvent(), ZeePlots::fillEleClassesPlots(), and ZIterativeAlgorithmWithFit::getWeight().

24 { return theElectron_; }
const reco::GsfElectron * theElectron_
Definition: CalibElectron.h:29

Member Data Documentation

const EcalRecHitCollection* calib::CalibElectron::theEEHits_
private

Definition at line 32 of file CalibElectron.h.

Referenced by getCalibModulesWeights(), and getEERecHits().

const reco::GsfElectron* calib::CalibElectron::theElectron_
private

Definition at line 29 of file CalibElectron.h.

Referenced by getCalibModulesWeights(), and getRecoElectron().

const EcalRecHitCollection* calib::CalibElectron::theHits_
private

Definition at line 31 of file CalibElectron.h.

Referenced by getCalibModulesWeights(), and getRecHits().