CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RxCalculator.cc
Go to the documentation of this file.
2 
4 
8 
12 
13 
14 using namespace edm;
15 using namespace reco;
16 
17 #define PI 3.141592653589793238462643383279502884197169399375105820974945
18 
19 
21 {
23  iEvent.getByLabel(hfLabel, hfhandle);
24  fHFRecHits_ = hfhandle.product();
25 
27  iEvent.getByLabel(hoLabel, hohandle);
28  fHORecHits_ = hohandle.product();
29 
31  iEvent.getByLabel(hbheLabel, hehbhandle);
32  fHBHERecHits_ = hehbhandle.product();
33 
34  ESHandle<CaloGeometry> geometryHandle;
35  iSetup.get<CaloGeometryRecord>().get(geometryHandle);
36  geometry_ = geometryHandle.product();
37 
38 }
39 
40 
41 double RxCalculator::getRx(const reco::SuperClusterRef cluster, double x, double threshold )
42 {
43  using namespace edm;
44  using namespace reco;
45 
46 
47  if(!fHBHERecHits_) {
48  LogError("RxCalculator") << "Error! Can't get HBHERecHits for event.";
49  return -100;
50  }
51 
52  double SClusterEta = cluster->eta();
53  double SClusterPhi = cluster->phi();
54  double TotalEt = 0;
55 
56  for(size_t index = 0; index < fHBHERecHits_->size(); index++) {
57  const HBHERecHit &rechit = (*fHBHERecHits_)[index];
58  const DetId &detid = rechit.id();
59  const GlobalPoint& hitpoint = geometry_->getPosition(detid);
60  double eta = hitpoint.eta();
61  double phi = hitpoint.phi();
62  double dEta = fabs(eta-SClusterEta);
63  double dPhi = fabs(phi-SClusterPhi);
64  while (dPhi>2*PI) dPhi-=2*PI;
65  if (dPhi>PI) dPhi=2*PI-dPhi;
66 
67  if (dPhi>PI) dPhi=2*PI-dPhi;
68 
69  double dR = sqrt(dEta * dEta + dPhi * dPhi);
70 
71  if (dR<x*0.1) {
72  double et = rechit.energy()/cosh(eta);
73  if (et<threshold) et=0;
74  TotalEt += et;
75  }
76  }
77 
78  return TotalEt;
79 }
80 
81 double RxCalculator::getROx(const reco::SuperClusterRef cluster, double x,double threshold)
82 {
83  using namespace edm;
84  using namespace reco;
85 
86  if(!fHORecHits_) {
87  LogError("RxCalculator") << "Error! Can't get HORecHits for event.";
88  return -100;
89  }
90 
91  double SClusterEta = cluster->eta();
92  double SClusterPhi = cluster->phi();
93  double TotalEt = 0;
94 
95  for(size_t index = 0; index < fHORecHits_->size(); index++) {
96  const HORecHit &rechit = (*fHORecHits_)[index];
97  const DetId &detid = rechit.id();
98  const GlobalPoint& hitpoint = geometry_->getPosition(detid);
99  double eta = hitpoint.eta();
100  double phi = hitpoint.phi();
101  double dEta = fabs(eta-SClusterEta);
102  double dPhi = fabs(phi-SClusterPhi);
103  while (dPhi>2*PI) dPhi-=2*PI;
104  if (dPhi>PI) dPhi=2*PI-dPhi;
105 
106  double dR = sqrt(dEta * dEta + dPhi * dPhi);
107  if (dR<x*0.1) {
108  double et = rechit.energy()/cosh(eta);
109  if (et<threshold) et=0;
110  TotalEt += et;
111  }
112  }
113  return TotalEt;
114 }
115 
116 double RxCalculator::getRFx(const reco::SuperClusterRef cluster, double x, double threshold)
117 {
118  using namespace edm;
119  using namespace reco;
120 
121  if(!fHFRecHits_) {
122  LogError("RxCalculator") << "Error! Can't get HFRecHits for event.";
123  return -100;
124  }
125 
126  double SClusterEta = cluster->eta();
127  double SClusterPhi = cluster->phi();
128  double TotalEt = 0;
129 
130  for(size_t index = 0; index < fHFRecHits_->size(); index++) {
131  const HFRecHit &rechit = (*fHFRecHits_)[index];
132  const DetId &detid = rechit.id();
133  const GlobalPoint& hitpoint = geometry_->getPosition(detid);
134  double eta = hitpoint.eta();
135  double phi = hitpoint.phi();
136  double dEta = fabs(eta-SClusterEta);
137  double dPhi = fabs(phi-SClusterPhi);
138  while (dPhi>2*PI) dPhi-=2*PI;
139  if (dPhi>PI) dPhi=2*PI-dPhi;
140 
141 
142  double dR = sqrt(dEta * dEta + dPhi * dPhi);
143  if (dR<x*0.1) {
144  double et = rechit.energy()/cosh(eta);
145  if (et<threshold) et=0;
146  TotalEt += et;
147  }
148  }
149 
150 
151 
152  return TotalEt;
153 }
154 
155 
156 double RxCalculator::getCRx(const reco::SuperClusterRef cluster, double x, double threshold )
157 {
158  using namespace edm;
159  using namespace reco;
160 
161 
162  if(!fHBHERecHits_) {
163  LogError("RxCalculator") << "Error! Can't get HBHERecHits for event.";
164  return -100;
165  }
166 
167  double SClusterEta = cluster->eta();
168  double SClusterPhi = cluster->phi();
169  double TotalEt = 0;
170 
171  for(size_t index = 0; index < fHBHERecHits_->size(); index++) {
172  const HBHERecHit &rechit = (*fHBHERecHits_)[index];
173  const DetId &detid = rechit.id();
174  const GlobalPoint& hitpoint = geometry_->getPosition(detid);
175  double eta = hitpoint.eta();
176  double phi = hitpoint.phi();
177  double dEta = fabs(eta-SClusterEta);
178  double dPhi = fabs(phi-SClusterPhi);
179  while (dPhi>2*PI) dPhi-=2*PI;
180  if (dPhi>PI) dPhi=2*PI-dPhi;
181 
182  if (dEta<x*0.1) {
183  double et = rechit.energy()/cosh(eta);
184  if (et<threshold) et=0;
185  TotalEt += et;
186  }
187  }
188 
189  double Rx = getRx(cluster,x,threshold);
190  double CRx = Rx - TotalEt / 40.0 * x;
191 
192  return CRx;
193 }
#define PI
HcalDetId id() const
get the id
Definition: HBHERecHit.h:21
double getCRx(const reco::SuperClusterRef clus, double i, double threshold)
Geom::Phi< T > phi() const
Definition: PV3DBase.h:63
T eta() const
int iEvent
Definition: GenABIO.cc:243
double dPhi(double phi1, double phi2)
Definition: JetUtil.h:30
float energy() const
Definition: CaloRecHit.h:19
T sqrt(T t)
Definition: SSEVec.h:28
double getRx(const reco::SuperClusterRef clus, double i, double threshold)
Definition: RxCalculator.cc:41
RxCalculator(const edm::Event &iEvent, const edm::EventSetup &iSetup, edm::InputTag hbheLabel, edm::InputTag hfLabel, edm::InputTag hoLabel)
Definition: RxCalculator.cc:20
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:355
Definition: DetId.h:20
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
T const * product() const
Definition: Handle.h:74
double getROx(const reco::SuperClusterRef clus, double i, double threshold)
Definition: RxCalculator.cc:81
T eta() const
Definition: PV3DBase.h:70
double getRFx(const reco::SuperClusterRef clus, double i, double threshold)
HcalDetId id() const
get the id
Definition: HFRecHit.h:21
Definition: DDAxes.h:10
HcalDetId id() const
get the id
Definition: HORecHit.h:21
Definition: DDAxes.h:10