CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MuonHOAcceptance.cc
Go to the documentation of this file.
2 
3 #include <iostream>
4 #include <algorithm>
5 #include <cmath>
6 
7 //#include "TTree.h"
8 #include "TMultiGraph.h"
9 #include "TGraph.h"
10 
12 
15 
18 
20 
21 std::vector<uint32_t> MuonHOAcceptance::deadIds;
22 std::vector<MuonHOAcceptance::deadRegion> MuonHOAcceptance::deadRegions;
23 std::vector<uint32_t> MuonHOAcceptance::SiPMIds;
24 std::vector<MuonHOAcceptance::deadRegion> MuonHOAcceptance::SiPMRegions;
25 bool MuonHOAcceptance::inited = false;
26 int const MuonHOAcceptance::etaBounds = 5;
27 double const MuonHOAcceptance::etaMin[etaBounds] =
28  //{-1.262, -0.861, -0.307, 0.341, 0.885};
29  {-1.2544, -0.8542, -0.3017, 0.3425, 0.8796};
30 double const MuonHOAcceptance::etaMax[etaBounds] =
31  //{-0.885, -0.341, 0.307, 0.861, 1.262};
32  {-0.8796, -0.3425, 0.3017, 0.8542, 1.2544};
33 double const MuonHOAcceptance::twopi = 2.*3.14159265358979323846;
34 int const MuonHOAcceptance::phiSectors = 12;
35 double const MuonHOAcceptance::phiMinR0[phiSectors] = {-0.16172,
36  0.3618786,
37  0.8854773,
38  1.409076116,
39  1.932674892,
40  2.456273667,
41  2.979872443,
42  3.503471219,
43  4.027069994,
44  4.55066877,
45  5.074267545,
46  5.597866321 };
47 double const MuonHOAcceptance::phiMaxR0[phiSectors] = { 0.317395374,
48  0.84099415,
49  1.364592925,
50  1.888191701,
51  2.411790477,
52  2.935389252,
53  3.458988028,
54  3.982586803,
55  4.506185579,
56  5.029784355,
57  5.55338313,
58  6.076981906 };
59 double const MuonHOAcceptance::phiMinR12[phiSectors] = {-0.166264081,
60  0.357334694,
61  0.88093347,
62  1.404532245,
63  1.928131021,
64  2.451729797,
65  2.975328572,
66  3.498927348,
67  4.022526123,
68  4.546124899,
69  5.069723674,
70  5.59332245 };
71 double const MuonHOAcceptance::phiMaxR12[phiSectors] = { 0.34398862,
72  0.867587396,
73  1.391186172,
74  1.914784947,
75  2.438383723,
76  2.961982498,
77  3.485581274,
78  4.00918005,
79  4.532778825,
80  5.056377601,
81  5.579976376,
82  6.103575152 };
83 
85  if (!inited) return false;
86  std::vector<uint32_t>::const_iterator found =
87  std::find(deadIds.begin(), deadIds.end(), id);
88  if ((found != deadIds.end()) && (*found == id)) return true;
89  else return false;
90 }
91 
93  if (!inited) return false;
94  std::vector<uint32_t>::const_iterator found =
95  std::find(SiPMIds.begin(), SiPMIds.end(), id);
96  if ((found != SiPMIds.end()) && (*found == id)) return true;
97  else return false;
98 }
99 
101  double delta_eta, double delta_phi)
102 {
103  for (int ieta = 0; ieta<etaBounds; ++ieta) {
104  if ( (eta > etaMin[ieta]+delta_eta) &&
105  (eta < etaMax[ieta]-delta_eta) ) {
106  for (int iphi = 0; iphi<phiSectors; ++iphi) {
107  double const * mins = ((ieta == 2) ? phiMinR0 : phiMinR12);
108  double const * maxes = ((ieta == 2) ? phiMaxR0 : phiMaxR12);
109  while (phi < mins[0])
110  phi += twopi;
111  while (phi > mins[0]+twopi)
112  phi -= twopi;
113  if ( ( phi > mins[iphi] + delta_phi ) &&
114  ( phi < maxes[iphi] - delta_phi ) ) {
115  return true;
116  }
117  }
118  return false;
119  }
120  }
121  return false;
122 }
123 
125  double delta_eta, double delta_phi) {
126  if (!inited) return true;
127  int ieta = int(eta/0.087) + ((eta>0) ? 1 : -1);
128  double const * mins = ((std::abs(ieta) > 4) ? phiMinR12 : phiMinR0);
129  while (phi < mins[0])
130  phi += twopi;
131  while (phi > mins[0]+twopi)
132  phi -= twopi;
133  std::vector<deadRegion>::const_iterator region;
134  for (region = deadRegions.begin(); region != deadRegions.end(); ++region) {
135  if ( (phi < region->phiMax + delta_phi) &&
136  (phi > region->phiMin - delta_phi) &&
137  (eta < region->etaMax + delta_eta) &&
138  (eta > region->etaMin - delta_eta) )
139  return false;
140  }
141  return true;
142 }
143 
144 bool MuonHOAcceptance::inSiPMGeom(double eta, double phi,
145  double delta_eta, double delta_phi) {
146  if (!inited) return false;
147  int ieta = int(eta/0.087) + ((eta>0) ? 1 : -1);
148  double const * mins = ((std::abs(ieta) > 4) ? phiMinR12 : phiMinR0);
149  while (phi < mins[0])
150  phi += twopi;
151  while (phi > mins[0]+twopi)
152  phi -= twopi;
153  std::vector<deadRegion>::const_iterator region;
154  for (region = SiPMRegions.begin(); region != SiPMRegions.end(); ++region) {
155  if ( (phi < region->phiMax - delta_phi) &&
156  (phi > region->phiMin + delta_phi) &&
157  (eta < region->etaMax - delta_eta) &&
158  (eta > region->etaMin + delta_eta) ) {
159  return true;
160  }
161  }
162  return false;
163 }
164 
166  deadIds.clear();
167 
169  eSetup.get<HcalChannelQualityRcd>().get(p);
170  HcalChannelQuality *myqual = new HcalChannelQuality(*p.product());
171 
173  eSetup.get<HcalSeverityLevelComputerRcd>().get(mycomputer);
174  const HcalSeverityLevelComputer *mySeverity = mycomputer.product();
175 
176  // TTree * deads = new TTree("deads", "deads");
177  // deads->ReadFile("HOdeadnessChannels.txt", "ieta/I:iphi/I:deadness/D");
178  int ieta, iphi;
179  // double deadness;
180  // deads->SetBranchAddress("ieta", &ieta);
181  // deads->SetBranchAddress("iphi", &iphi);
182  // deads->SetBranchAddress("deadness", &deadness);
183  // deads->Print();
184  //std::cout << "ieta\tiphi\n";
185  for (ieta=-15; ieta <= 15; ieta++) {
186  if (ieta != 0) {
187  for (iphi = 1; iphi <= 72; iphi++) {
188  // for (int i=0; i<deads->GetEntries(); ++i) {
189  // deads->GetEntry(i);
190  // if (deadness > 0.4) {
191  HcalDetId did(HcalOuter,ieta,iphi,4);
192  const HcalChannelStatus *mystatus = myqual->getValues(did.rawId());
193  if (mySeverity->dropChannel(mystatus->getValue())) {
194  deadIds.push_back(did.rawId());
195  // std::cout << did.ieta() << '\t' << did.iphi() << '\n';
196  }
197  //HO +1 RBX 10
198  if ( (ieta>=5) && (ieta<=10) && (iphi >= 47) && (iphi <= 58) ) {
199  SiPMIds.push_back(did.rawId());
200  }
201  //HO +2 RBX 12
202  if ( (ieta>=11) && (ieta<=15) && (iphi >= 59) && (iphi <= 70) ) {
203  SiPMIds.push_back(did.rawId());
204  }
205  }
206  }
207  }
208  std::sort(deadIds.begin(), deadIds.end());
209  std::sort(SiPMIds.begin(), SiPMIds.end());
210  // std::cout << "SiPMIds: " << SiPMIds.size() << '\n';
211  // delete deads;
212  buildDeadAreas();
213  buildSiPMAreas();
214  inited = true;
215  delete myqual;
216 }
217 
219  std::vector<uint32_t>::iterator did;
220  std::list<deadIdRegion> didregions;
221  for (did = deadIds.begin(); did != deadIds.end(); ++did) {
222  HcalDetId tmpId(*did);
223  didregions.push_back( deadIdRegion( tmpId.ieta(), tmpId.ieta(),
224  tmpId.iphi(), tmpId.iphi() ) );
225  }
226  // std::cout << "dead regions: " << didregions.size() << '\n';
227 
228  mergeRegionLists(didregions);
229  convertRegions(didregions, deadRegions);
230 }
231 
233  std::vector<uint32_t>::iterator sid;
234  std::list<deadIdRegion> idregions;
235 
236  for (sid = SiPMIds.begin(); sid != SiPMIds.end(); ++sid) {
237  HcalDetId tmpId(*sid);
238  idregions.push_back( deadIdRegion( tmpId.ieta(), tmpId.ieta(),
239  tmpId.iphi(), tmpId.iphi() ) );
240  }
241 
242  mergeRegionLists(idregions);
243  convertRegions(idregions,SiPMRegions);
244 }
245 
246 void MuonHOAcceptance::mergeRegionLists (std::list<deadIdRegion>& didregions) {
247  std::list<deadIdRegion>::iterator curr;
248  std::list<deadIdRegion> list2;
249  unsigned int startSize;
250  do {
251  startSize = didregions.size();
252 
253  // std::cout << "regions: " << startSize << '\n';
254  //merge in phi
255  curr = didregions.begin();
256  while (curr != didregions.end()) {
257  deadIdRegion merger(*curr);
258  curr = didregions.erase(curr);
259  while (curr != didregions.end()) {
260  if ( (merger.sameEta(*curr)) && (merger.adjacentPhi(*curr)) ) {
261  merger.merge(*curr);
262  curr = didregions.erase(curr);
263  } else ++curr;
264  }
265  list2.push_back(merger);
266  curr = didregions.begin();
267  }
268 
269  //merge in eta
270  curr = list2.begin();
271  while (curr != list2.end()) {
272  deadIdRegion merger(*curr);
273  curr = list2.erase(curr);
274  while (curr != list2.end()) {
275  if ( (merger.samePhi(*curr)) && (merger.adjacentEta(*curr)) ) {
276  merger.merge(*curr);
277  curr = list2.erase(curr);
278  } else ++curr;
279  }
280  didregions.push_back(merger);
281  curr = list2.begin();
282  }
283  } while (startSize > didregions.size());
284 }
285 
286 void MuonHOAcceptance::convertRegions(std::list<deadIdRegion> const& idregions,
287  std::vector<deadRegion>& regions) {
288  double e1, e2;
289  double eMin,eMax,pMin,pMax;
290  static double const etaStep = 0.087;
291  static double const phiStep = twopi/72.;
292  double const offset = 2.;
293  double const * mins;
294  double const * maxes;
295  std::list<deadIdRegion>::const_iterator curr;
296  double zero;
297  for (curr = idregions.begin(); curr != idregions.end(); ++curr) {
298  // std::cout << "region boundaries: ieta,iphi\n"
299  // << " min: " << curr->etaMin << ',' << curr->phiMin << '\n'
300  // << " max: " << curr->etaMax << ',' << curr->phiMax << '\n';
301  if (curr->etaMin == -4) {
302  eMin = etaMax[1];
303  } else if (curr->etaMin == 5) {
304  eMin = etaMin[3];
305  } else {
306  e1 = (std::abs(curr->etaMin)-1)*etaStep*
307  (-(curr->etaMin<0) + (curr->etaMin>0));
308  e2 = std::abs(curr->etaMin)*etaStep*
309  (-(curr->etaMin<0) + (curr->etaMin>0));
310  eMin = std::min(e1,e2);
311  }
312  if (curr->etaMax == 4) {
313  eMax = etaMin[3];
314  } else if (curr->etaMax == -5) {
315  eMax = etaMax[1];
316  } else {
317  e1 = (std::abs(curr->etaMax)-1)*etaStep*
318  (-(curr->etaMax<0) + (curr->etaMax>0));
319  e2 = std::abs(curr->etaMax)*etaStep*
320  (-(curr->etaMax<0) + (curr->etaMax>0));
321  eMax = std::max(e1,e2);
322  }
323  mins = ((std::abs(curr->etaMin)>4) ? phiMinR12 : phiMinR0);
324  maxes = ((std::abs(curr->etaMin)>4) ? phiMaxR12 : phiMaxR0);
325  zero = (mins[0] + maxes[phiSectors-1] - twopi)/2.;
326  pMin = (curr->phiMin-1)*phiStep + zero + phiStep*offset;
327  pMax = curr->phiMax*phiStep + zero + phiStep*offset;
328  while (pMax < mins[0])
329  pMax += twopi;
330  while (pMax > mins[0]+twopi)
331  pMax -= twopi;
332  while (pMin < mins[0])
333  pMin += twopi;
334  while (pMin > mins[0]+twopi)
335  pMin -= twopi;
336 
337  regions.push_back( deadRegion(eMin, eMax, pMin, pMax) );
338  // std::cout << " : eta,phi\n"
339  // << " min: " << tmp.etaMin << ',' << tmp.phiMin << '\n'
340  // << " max: " << tmp.etaMax << ',' << tmp.phiMax << '\n';
341  }
342 }
343 
344 TMultiGraph * MuonHOAcceptance::graphRegions(std::vector<deadRegion> const& regions) {
345  TMultiGraph * bounds = new TMultiGraph("bounds", "bounds");
346  std::vector<deadRegion>::const_iterator region;
347  TGraph * gr;
348  for (region = regions.begin(); region != regions.end(); ++region) {
349  // std::cout << "region eta range:(" << region->etaMin << ',' << region->etaMax
350  // << ") phi range: (" << region->phiMin << ',' << region->phiMax << ")\n";
351  double pMin = region->phiMin;
352  while (pMin > twopi/2.) pMin -= twopi;
353  double pMax = region->phiMax;
354  while (pMax > twopi/2.) pMax -= twopi;
355 
356  if (pMin < pMax) {
357  gr = new TGraph(5);
358  gr->SetLineColor(kRed);
359  gr->SetPoint(0, region->etaMin, pMin);
360  gr->SetPoint(1, region->etaMin, pMax);
361  gr->SetPoint(2, region->etaMax, pMax);
362  gr->SetPoint(3, region->etaMax, pMin);
363  gr->SetPoint(4, region->etaMin, pMin);
364  bounds->Add(gr, "l");
365  } else {
366  gr = new TGraph(5);
367  gr->SetLineColor(kRed);
368  gr->SetPoint(0, region->etaMin, pMin);
369  gr->SetPoint(1, region->etaMin, twopi/2.);
370  gr->SetPoint(2, region->etaMax, twopi/2.);
371  gr->SetPoint(3, region->etaMax, pMin);
372  gr->SetPoint(4, region->etaMin, pMin);
373  bounds->Add(gr, "l");
374  gr = new TGraph(5);
375  gr->SetLineColor(kRed);
376  gr->SetPoint(0, region->etaMin, -twopi/2.);
377  gr->SetPoint(1, region->etaMin, pMax);
378  gr->SetPoint(2, region->etaMax, pMax);
379  gr->SetPoint(3, region->etaMax, -twopi/2.);
380  gr->SetPoint(4, region->etaMin, -twopi/2.);
381  bounds->Add(gr, "l");
382  }
383  }
384  return bounds;
385 }
386 
388  etaMin = std::min(etaMin, other.etaMin);
389  etaMax = std::max(etaMax, other.etaMax);
390  phiMin = std::min(phiMin, other.phiMin);
391  phiMax = std::max(phiMax, other.phiMax);
392 }
bool adjacentPhi(deadIdRegion const &other)
static double const twopi
static bool inSiPMGeom(double eta, double phi, double delta_eta=0., double delta_phi=0.)
double delta_eta(double eta1, double eta2)
Definition: AnglesUtil.h:98
bool sameEta(deadIdRegion const &other)
static TMultiGraph * graphRegions(std::vector< deadRegion > const &regions)
static void mergeRegionLists(std::list< deadIdRegion > &didregions)
static std::vector< uint32_t > deadIds
static double const phiMaxR0[]
static double const phiMinR12[]
const Item * getValues(DetId fId, bool throwOnFail=true) const
T eta() const
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
static void initIds(edm::EventSetup const &eSetup)
static int const etaBounds
static double const phiMinR0[]
bool samePhi(deadIdRegion const &other)
const T & max(const T &a, const T &b)
static void buildDeadAreas()
int ieta() const
get the cell ieta
Definition: HcalDetId.h:36
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool dropChannel(const uint32_t &mystatus) const
unsigned int offset(bool)
static double const etaMax[]
static bool inNotDeadGeom(double eta, double phi, double delta_eta=0., double delta_phi=0.)
double delta_phi(double ph11, double phi2)
Definition: AnglesUtil.h:91
bool adjacentEta(deadIdRegion const &other)
int iphi() const
get the cell iphi
Definition: HcalDetId.h:38
Float e1
Definition: deltaR.h:22
static void buildSiPMAreas()
static std::vector< deadRegion > deadRegions
static double const etaMin[]
static double const phiMaxR12[]
const T & get() const
Definition: EventSetup.h:55
static bool isChannelDead(uint32_t id)
T const * product() const
Definition: ESHandle.h:62
static bool inited
Float e2
Definition: deltaR.h:23
static bool inGeomAccept(double eta, double phi, double delta_eta=0., double delta_phi=0.)
static int const phiSectors
static std::vector< uint32_t > SiPMIds
static void convertRegions(std::list< deadIdRegion > const &idregions, std::vector< deadRegion > &regions)
void merge(deadIdRegion const &other)
uint32_t getValue() const
static std::vector< deadRegion > SiPMRegions
static bool isChannelSiPM(uint32_t id)
Definition: DDAxes.h:10