CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_9/src/CondCore/HcalPlugins/interface/HcalObjRepresent.h

Go to the documentation of this file.
00001 #ifndef HcalObjRepresent_h
00002 #define HcalObjRepresent_h
00003 
00004 #include "CondFormats/HcalObjects/interface/HcalGains.h"
00005 
00006 //#include "CondCore/Utilities/interface/PayLoadInspector.h"
00007 //#include "CondCore/Utilities/interface/InspectorPythonWrapper.h"
00008 
00009 #include <string>
00010 #include <fstream>
00011 #include <sstream>
00012 #include <vector>
00013 
00014 #include "TH1F.h"
00015 #include "TH2F.h"
00016 #include "DataFormats/DetId/interface/DetId.h"
00017 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
00018 
00019 #include "TROOT.h"
00020 #include "TCanvas.h"
00021 #include "TStyle.h"
00022 #include "TColor.h"
00023 #include "TLine.h"
00024 
00025 //functions for correct representation of data in summary and plot
00026 namespace HcalObjRepresent{
00027         inline std::string IntToBinary(unsigned int number) {
00028                 std::stringstream ss;
00029                 unsigned int mask = 1<<31;
00030                 for (unsigned short int i = 0; i < 32; ++i){
00031                         //if (!(i % 4))
00032                         //      ss << "_";
00033                         if (mask & number)
00034                                 ss << "1";
00035                         else 
00036                                 ss << "0";
00037                         mask = mask >> 1;
00038                 }
00039                 return ss.str();
00040         }
00041 
00042 
00043         const bool isBitSet(unsigned int bitnumber, unsigned int status)
00044         { 
00045                 unsigned int statadd = 0x1<<(bitnumber);
00046                 return (status&statadd)?(true):(false);
00047         }
00048 
00049 
00050         std::string getBitsSummary(uint32_t bits, std::string  statusBitArray[], short unsigned int bitMap[]  ){
00051                 std::stringstream ss;
00052                 for (unsigned int i = 0; i < 9; ++i){
00053                         if (isBitSet(bitMap[i], bits)){
00054                                 ss << "[" <<bitMap[i]<< "]" << statusBitArray[bitMap[i]] << "; ";
00055                         }
00056                 }
00057                 ss << std::endl;
00058                 return ss.str();
00059         }
00060 
00061 
00062         //functions for making plot:
00063         void setBinLabels(std::vector<TH2F> &depth)
00064         {
00065                 // Set labels for all depth histograms
00066                 for (unsigned int i=0;i<depth.size();++i)
00067                 {
00068                         depth[i].SetXTitle("i#eta");
00069                         depth[i].SetYTitle("i#phi");
00070                 }
00071 
00072                 std::stringstream label;
00073 
00074                 // set label on every other bin
00075                 for (int i=-41;i<=-29;i=i+2)
00076                 {
00077                         label<<i;
00078                         depth[0].GetXaxis()->SetBinLabel(i+42,label.str().c_str());
00079                         depth[1].GetXaxis()->SetBinLabel(i+42,label.str().c_str());
00080                         label.str("");
00081                 }
00082                 depth[0].GetXaxis()->SetBinLabel(14,"-29HE");
00083                 depth[1].GetXaxis()->SetBinLabel(14,"-29HE");
00084 
00085                 // offset by one for HE
00086                 for (int i=-27;i<=27;i=i+2)
00087                 {
00088                         label<<i;
00089                         depth[0].GetXaxis()->SetBinLabel(i+43,label.str().c_str());
00090                         label.str("");
00091                 }
00092                 depth[0].GetXaxis()->SetBinLabel(72,"29HE");
00093                 for (int i=29;i<=41;i=i+2)
00094                 {
00095                         label<<i;
00096                         depth[0].GetXaxis()->SetBinLabel(i+44,label.str().c_str());
00097                         label.str("");
00098                 }
00099                 for (int i=16;i<=28;i=i+2)
00100                 {
00101                         label<<i-43;
00102                         depth[1].GetXaxis()->SetBinLabel(i,label.str().c_str());
00103                         label.str("");
00104                 }
00105                 depth[1].GetXaxis()->SetBinLabel(29,"NULL");
00106                 for (int i=15;i<=27;i=i+2)
00107                 {
00108                         label<<i;
00109                         depth[1].GetXaxis()->SetBinLabel(i+15,label.str().c_str());
00110                         label.str("");
00111                 }
00112 
00113                 depth[1].GetXaxis()->SetBinLabel(44,"29HE");
00114                 for (int i=29;i<=41;i=i+2)
00115                 {
00116                         label<<i;
00117                         depth[1].GetXaxis()->SetBinLabel(i+16,label.str().c_str());
00118                         label.str("");
00119                 }
00120 
00121                 // HE depth 3 labels;
00122                 depth[2].GetXaxis()->SetBinLabel(1,"-28");
00123                 depth[2].GetXaxis()->SetBinLabel(2,"-27");
00124                 depth[2].GetXaxis()->SetBinLabel(3,"Null");
00125                 depth[2].GetXaxis()->SetBinLabel(4,"-16");
00126                 depth[2].GetXaxis()->SetBinLabel(5,"Null");
00127                 depth[2].GetXaxis()->SetBinLabel(6,"16");
00128                 depth[2].GetXaxis()->SetBinLabel(7,"Null");
00129                 depth[2].GetXaxis()->SetBinLabel(8,"27");
00130                 depth[2].GetXaxis()->SetBinLabel(9,"28");
00131         }
00132         // Now define functions that can be used in conjunction with EtaPhi histograms
00133 
00134         // This arrays the eta binning for depth 2 histograms (with a gap between -15 -> +15)
00135         const int binmapd2[]={-42,-41,-40,-39,-38,-37,-36,-35,-34,-33,-32,-31,-30,
00136                 -29,-28,-27,-26,-25,-24,-23,-22,-21,-20,-19,-18,-17,
00137                 -16,-15,-9999, 15,16,17,18,19,20,21,22,23,24,25,26,27,28,29,
00138                 30,31,32,33,34,35,36,37,38,39,40,41,42};
00139 
00140         // This stores eta binning in depth 3 (where HE is only present at a few ieta values)
00141 
00142         const int binmapd3[]={-28,-27,-9999,-16,-9999,16,-9999,27,28};
00143 
00144         inline int CalcEtaBin(int subdet, int ieta, int depth)
00145         {
00146                 // This takes the eta value from a subdetector and return an eta counter value as used by eta-phi array
00147                 // (ieta=-41 corresponds to bin 0, +41 to bin 85 -- there are two offsets to deal with the overlap at |ieta|=29).
00148                 // For HO, ieta = -15 corresponds to bin 0, and ieta=15 is bin 30
00149                 // For HE depth 3, things are more complicated, but feeding the ieta value will give back the corresponding counter eta value
00150 
00151                 // The CalcEtaBin value is the value as used within our array counters, and thus starts at 0.
00152                 // If you are using it with getBinContent or setBinContent, you will need to add +1 to the result of this function
00153 
00154                 int etabin=-9999; // default invalid value
00155 
00156                 if (depth==1)
00157                 {
00158                         // Depth 1 is fairly straightforward -- just shift HF-, HF+ by -/+1
00159                         etabin=ieta+42;
00160                         if (subdet==HcalForward)
00161                         {
00162                                 ieta < 0 ? etabin-- : etabin++;
00163                         }
00164                 }
00165 
00166                 else if (depth==2)
00167                 {
00168                         // Depth 2 is more complicated, given that there are no cells in the range |ieta|<15
00169                         if (ieta<-14)
00170                         {
00171                                 etabin=ieta+42;
00172                                 if (subdet==HcalForward) etabin--;
00173                         }
00174                         else if (ieta>14)
00175                         {
00176                                 etabin=ieta+14;
00177                                 if (subdet==HcalForward) etabin++;
00178                         }
00179 
00180                 }
00181                 // HO is also straightforward; a simple offset to the ieta value is applied
00182                 else if (subdet==HcalOuter && abs(ieta)<16)
00183                         etabin=ieta+15;
00184                 else if (subdet==HcalEndcap)
00185                 {
00186                         // HE depth 3 has spotty coverage; hard-code the bin response
00187                         if (depth==3)
00188                         {
00189                                 if (ieta==-28) etabin=0;
00190                                 else if (ieta==-27) etabin=1;
00191                                 else if (ieta==-16) etabin=3;
00192                                 else if (ieta==16)  etabin=5;
00193                                 else if (ieta==27)  etabin=7;
00194                                 else if (ieta==28)  etabin=8;
00195                         }
00196                 }
00197                 return etabin;
00198         }
00199 
00200         inline int CalcIeta(int subdet, int eta, int depth)
00201         {
00202                 // This function returns the 'true' ieta value given subdet, eta, and depth
00203                 // Here 'eta' is the index from our arrays (it starts at 0);
00204                 // remember that histogram bins start with bin 1, so there's an offset of 1
00205                 // to consider if using getBinContent(eta,phi)
00206 
00207                 // eta runs from 0...X  (X depends on depth)
00208                 int ieta=-9999; // default value is nonsensical
00209                 if (subdet==HcalBarrel)
00210                 {
00211                         if (depth==1) 
00212                         {
00213                                 ieta=eta-42;
00214                                 if (ieta==0) return -9999;
00215                                 return ieta;
00216                         }
00217                         else if (depth==2)
00218                         {
00219                                 ieta=binmapd2[eta];
00220                                 if (ieta==0) return -9999;
00221                                 if (ieta==17 || ieta == -17) 
00222                                         return -9999; // no depth 2 cells at |ieta| = 17
00223                                 return ieta;
00224                         }
00225                         else
00226                                 return -9999; // non-physical value
00227                 }
00228                 else if (subdet==HcalForward)
00229                 {
00230                         if (depth==1)
00231                         {
00232                                 ieta=eta-42;
00233                                 if (eta<13) ieta++;
00234                                 else if (eta>71) ieta--;
00235                                 else return -9999; // if outside forward range, return dummy
00236                                 return ieta;
00237                         }
00238                         else if (depth==2)
00239                         {
00240                                 ieta=binmapd2[eta]; // special map for depth 2
00241                                 if (ieta<=-30) ieta++;
00242                                 else if (ieta>=30) ieta--;
00243                                 else return -9999;
00244                                 return ieta;
00245                         }
00246                         else return -9999;
00247                 }
00248 
00249                 else if (subdet==HcalEndcap)
00250                 {
00251                         if (depth==1) 
00252                                 ieta=eta-42;
00253                         else if (depth==2) 
00254                         {
00255                                 ieta=binmapd2[eta];
00256                                 if (abs(ieta)>29 || abs(ieta)<18) return -9999; // outside HE
00257                                 if (ieta==0) return -9999;
00258                                 return ieta;
00259                         }
00260                         else if (depth==3)
00261                         {
00262                                 if (eta<0 || eta>8) return -9999;
00263                                 else
00264                                         ieta=binmapd3[eta]; // special map for depth 3
00265                                 if (ieta==0) return -9999;
00266                                 return ieta;
00267                         }
00268                         else return -9999;
00269                 } // HcalEndcap
00270                 else if ( subdet==HcalOuter)
00271                 {
00272                         if (depth!=4)
00273                                 return -9999;
00274                         else
00275                         {
00276                                 ieta= eta-15;  // bin 0 is ieta=-15, all bins increment normally from there
00277                                 if (abs(ieta)>15) return -9999;
00278                                 if (ieta==0) return -9999;
00279                                 return ieta;
00280                         }
00281                 } // HcalOuter
00282                 if (ieta==0) return -9999;
00283                 return ieta;
00284         }
00285 
00286         inline int CalcIeta(int eta, int depth)
00287         {
00288                 // This version of CalcIeta does the same as the function above,
00289                 // but does not require that 'subdet' be specified.
00290 
00291                 // returns ieta value give an eta counter.
00292                 // eta runs from 0...X  (X depends on depth)
00293                 int ieta=-9999;
00294                 if (eta<0) return ieta;
00295                 if (depth==1)
00296                 {
00297                         ieta=eta-42; // default shift: bin 0 corresponds to a histogram ieta of -42 (which is offset by 1 from true HF value of -41)
00298                         if (eta<13) ieta++;
00299                         else if (eta>71) ieta--;
00300                         if (ieta==0) ieta=-9999;
00301                         return ieta;
00302                 }
00303                 else if (depth==2)
00304                 {
00305                         if (eta>57) return -9999;
00306                         else
00307                         {
00308                                 ieta=binmapd2[eta];
00309                                 if (ieta==-9999) return ieta;
00310                                 if (ieta==0) return -9999;
00311                                 if (ieta==17 || ieta == -17) return -9999; // no depth 2 cells at |ieta| = 17
00312                                 else if (ieta<=-30) ieta++;
00313                                 else if (ieta>=30) ieta--;
00314                                 return ieta;
00315                         }
00316                 }
00317                 else if (depth==3)
00318                 {
00319                         if (eta>8) return -9999;
00320                         else
00321                                 ieta=binmapd3[eta];
00322                         if (ieta==0) return -9999;
00323                         return ieta;
00324                 }
00325                 else if (depth==4)
00326                 {
00327                         ieta= eta-15;  // bin 0 is ieta=-15, all bins increment normally from there
00328                         if (abs(ieta)>15) return -9999;
00329                         if (ieta==0) return -9999;
00330                         return ieta;
00331                 }
00332                 return ieta; // avoids compilation warning
00333         }
00334 
00335 
00336         // Functions to check whether a given (eta,depth) value is valid for a given subdetector
00337 
00338         inline std::vector<std::string> HcalEtaPhiHistNames()
00339         {
00340                 std::vector<std::string> name;
00341                 name.push_back("HB HE HF Depth 1 ");
00342                 name.push_back("HB HE HF Depth 2 ");
00343                 name.push_back("HE Depth 3 ");
00344                 name.push_back("HO Depth 4 ");
00345                 return name;
00346         }
00347 
00348 
00349         inline bool isHB(int etabin, int depth)
00350         {
00351                 if (depth>2) return false;
00352                 else if (depth<1) return false;
00353                 else
00354                 {
00355                         int ieta=CalcIeta(etabin,depth);
00356                         if (ieta==-9999) return false;
00357                         if (depth==1)
00358                         {
00359                                 if (abs(ieta)<=16 ) return true;
00360                                 else return false;
00361                         }
00362                         else if (depth==2)
00363                         {
00364                                 if (abs(ieta)==15 || abs(ieta)==16) return true;
00365                                 else return false;
00366                         }
00367                 }
00368                 return false;
00369         }
00370 
00371         inline bool isHE(int etabin, int depth)
00372         {
00373                 if (depth>3) return false;
00374                 else if (depth<1) return false;
00375                 else
00376                 {
00377                         int ieta=CalcIeta(etabin,depth);
00378                         if (ieta==-9999) return false;
00379                         if (depth==1)
00380                         {
00381                                 if (abs(ieta)>=17 && abs(ieta)<=28 ) return true;
00382                                 if (ieta==-29 && etabin==13) return true; // HE -29
00383                                 if (ieta==29 && etabin == 71) return true; // HE +29
00384                         }
00385                         else if (depth==2)
00386                         {
00387                                 if (abs(ieta)>=17 && abs(ieta)<=28 ) return true;
00388                                 if (ieta==-29 && etabin==13) return true; // HE -29
00389                                 if (ieta==29 && etabin == 43) return true; // HE +29
00390                         }
00391                         else if (depth==3)
00392                                 return true;
00393                 }
00394                 return false;
00395         }
00396 
00397         inline bool isHF(int etabin, int depth)
00398         {
00399                 if (depth>2) return false;
00400                 else if (depth<1) return false;
00401                 else
00402                 {
00403                         int ieta=CalcIeta(etabin,depth);
00404                         if (ieta==-9999) return false;
00405                         if (depth==1)
00406                         {
00407                                 if (ieta==-29 && etabin==13) return false; // HE -29
00408                                 else if (ieta==29 && etabin == 71) return false; // HE +29
00409                                 else if (abs(ieta)>=29 ) return true;
00410                         }
00411                         else if (depth==2)
00412                         {
00413                                 if (ieta==-29 && etabin==13) return false; // HE -29
00414                                 else if (ieta==29 && etabin==43) return false; // HE +29
00415                                 else if (abs(ieta)>=29 ) return true;
00416                         }
00417                 }
00418                 return false;
00419         }
00420 
00421         inline bool isHO(int etabin, int depth)
00422         {
00423                 if (depth!=4) return false;
00424                 int ieta=CalcIeta(etabin,depth);
00425                 if (ieta!=-9999) return true;
00426                 return false;
00427         }
00428 
00429         // Checks whether HO region contains SiPM
00430 
00431         inline bool isSiPM(int ieta, int iphi, int depth)
00432         {
00433                 if (depth!=4) return false;
00434                 // HOP1
00435                 if (ieta>=5 && ieta <=10 && iphi>=47 && iphi<=58) return true;  
00436                 // HOP2
00437                 if (ieta>=11 && ieta<=15 && iphi>=59 && iphi<=70) return true;
00438                 return false;
00439         }  // bool isSiPM
00440 
00441 
00442         // Checks whether (subdet, ieta, iphi, depth) value is a valid Hcal cell
00443 
00444         inline bool validDetId(HcalSubdetector sd, int ies, int ip, int dp)
00445         {
00446                 // inputs are (subdetector, ieta, iphi, depth)
00447                 // stolen from latest version of DataFormats/HcalDetId/src/HcalDetId.cc (not yet available in CMSSW_2_1_9)
00448 
00449                 const int ie ( abs( ies ) ) ;
00450 
00451                 return ( ( ip >=  1         ) &&
00452                         ( ip <= 72         ) &&
00453                         ( dp >=  1         ) &&
00454                         ( ie >=  1         ) &&
00455                         ( ( ( sd == HcalBarrel ) &&
00456                         ( ( ( ie <= 14         ) &&
00457                         ( dp ==  1         )    ) ||
00458                         ( ( ( ie == 15 ) || ( ie == 16 ) ) && 
00459                         ( dp <= 2          )                ) ) ) ||
00460                         (  ( sd == HcalEndcap ) &&
00461                         ( ( ( ie == 16 ) &&
00462                         ( dp ==  3 )          ) ||
00463                         ( ( ie == 17 ) &&
00464                         ( dp ==  1 )          ) ||
00465                         ( ( ie >= 18 ) &&
00466                         ( ie <= 20 ) &&
00467                         ( dp <=  2 )          ) ||
00468                         ( ( ie >= 21 ) &&
00469                         ( ie <= 26 ) &&
00470                         ( dp <=  2 ) &&
00471                         ( ip%2 == 1 )         ) ||
00472                         ( ( ie >= 27 ) &&
00473                         ( ie <= 28 ) &&
00474                         ( dp <=  3 ) &&
00475                         ( ip%2 == 1 )         ) ||
00476                         ( ( ie == 29 ) &&
00477                         ( dp <=  2 ) &&
00478                         ( ip%2 == 1 )         )          )      ) ||
00479                         (  ( sd == HcalOuter ) &&
00480                         ( ie <= 15 ) &&
00481                         ( dp ==  4 )           ) ||
00482                         (  ( sd == HcalForward ) &&
00483                         ( dp <=  2 )          &&
00484                         ( ( ( ie >= 29 ) &&
00485                         ( ie <= 39 ) &&
00486                         ( ip%2 == 1 )    ) ||
00487                         ( ( ie >= 40 ) &&
00488                         ( ie <= 41 ) &&
00489                         ( ip%4 == 3 )         )  ) ) ) ) ;
00490 
00491 
00492 
00493         } // bool validDetId(HcalSubdetector sd, int ies, int ip, int dp)
00494 
00495 
00496         // Sets eta, phi labels for 'summary' eta-phi plots (identical to Depth 1 Eta-Phi labelling)
00497 
00498         inline void SetEtaPhiLabels(TH2F &h)
00499         {
00500                 std::stringstream label;
00501                 for (int i=-41;i<=-29;i=i+2)
00502                 {
00503                         label<<i;
00504                         h.GetXaxis()->SetBinLabel(i+42,label.str().c_str());
00505                         label.str("");
00506                 }
00507                 h.GetXaxis()->SetBinLabel(14,"-29HE");
00508 
00509                 // offset by one for HE
00510                 for (int i=-27;i<=27;i=i+2)
00511                 {
00512                         label<<i;
00513                         h.GetXaxis()->SetBinLabel(i+43,label.str().c_str());
00514                         label.str("");
00515                 }
00516                 h.GetXaxis()->SetBinLabel(72,"29HE");
00517                 for (int i=29;i<=41;i=i+2)
00518                 {
00519                         label<<i;
00520                         h.GetXaxis()->SetBinLabel(i+44,label.str().c_str());
00521                         label.str("");
00522                 }
00523                 return;
00524         }
00525 
00526 
00527         // Fill Unphysical bins in histograms
00528         inline void FillUnphysicalHEHFBins(std::vector<TH2F> &hh)
00529         {
00530                 int ieta=0;
00531                 int iphi=0;
00532                 // First 2 depths have 5-10-20 degree corrections
00533                 for (unsigned int d=0;d<3;++d)
00534                 {
00535                         //BUG CAN BE HERE:
00536                         //if (hh[d] != 0) continue;
00537 
00538                         for (int eta=0;eta<hh[d].GetNbinsX();++eta)
00539                         {
00540                                 ieta=CalcIeta(eta,d+1);
00541                                 if (ieta==-9999 || abs(ieta)<21) continue;
00542                                 for (int phi=0;phi <hh[d].GetNbinsY();++phi)
00543                                 {
00544                                         iphi=phi+1;
00545                                         if (iphi%2==1 && abs(ieta)<40 && iphi<73)
00546                                         {
00547                                                 hh[d].SetBinContent(eta+1,iphi+1,hh[d].GetBinContent(eta+1,iphi));
00548                                         }
00549                                         // last two eta strips span 20 degrees in phi
00550                                         // Fill the phi cell above iphi, and the 2 below it
00551                                         else  if (abs(ieta)>39 && iphi%4==3 && iphi<73)
00552                                         {
00553                                                 //ieta=40, iphi=3 covers iphi 3,4,5,6
00554                                                 hh[d].SetBinContent(eta+1,(iphi)%72+1, hh[d].GetBinContent(eta+1,iphi));
00555                                                 hh[d].SetBinContent(eta+1,(iphi+1)%72+1, hh[d].GetBinContent(eta+1,iphi));
00556                                                 hh[d].SetBinContent(eta+1,(iphi+2)%72+1, hh[d].GetBinContent(eta+1,iphi));
00557                                         }
00558                                 } // for (int phi...)
00559                         } // for (int eta...)
00560                 } // for (int d=0;...)
00561                 // no corrections needed for HO (depth 4)
00562                 return;
00563         } // FillUnphysicalHEHFBins(MonitorElement* hh)
00564 
00565 
00566         //Fill unphysical bins for single ME
00567         inline void FillUnphysicalHEHFBins(TH2F &hh)
00568         {
00569                 // Fills unphysical HE/HF bins for Summary Histogram
00570                 // Summary Histogram is binned with the same binning as the Depth 1 EtaPhiHists
00571 
00572                 //CAN BE BUG HERE:
00573                 //if (hh==0) return;
00574 
00575                 int ieta=0;
00576                 int iphi=0;
00577                 int etabins = hh.GetNbinsX();
00578                 int phibins = hh.GetNbinsY();
00579                 float binval=0;
00580                 for (int eta=0;eta<etabins;++eta) // loop over eta bins
00581                 {
00582                         ieta=CalcIeta(eta,1);
00583                         if (ieta==-9999 || abs(ieta)<21) continue;  // ignore etas that don't exist, or that have 5 degree phi binning
00584 
00585                         for (int phi=0;phi<phibins;++phi)
00586                         {
00587                                 iphi=phi+1;
00588                                 if (iphi%2==1 && abs(ieta)<40 && iphi<73) // 10 degree phi binning condition
00589                                 {
00590                                         binval=hh.GetBinContent(eta+1,iphi);
00591                                         hh.SetBinContent(eta+1,iphi+1,binval);
00592                                 } // if (iphi%2==1...) 
00593                                 else if (abs(ieta)>39 && iphi%4==3 && iphi<73) // 20 degree phi binning condition
00594                                 {
00595                                         // Set last two eta strips where each cell spans 20 degrees in phi
00596                                         // Set next phi cell above iphi, and 2 cells below the actual cell 
00597                                         hh.SetBinContent(eta+1, (iphi)%72+1, hh.GetBinContent(eta+1,iphi));
00598                                         hh.SetBinContent(eta+1, (iphi+1)%72+1, hh.GetBinContent(eta+1,iphi));
00599                                         hh.SetBinContent(eta+1, (iphi+2)%72+1, hh.GetBinContent(eta+1,iphi));
00600                                 } // else if (abs(ieta)>39 ...)
00601                         } // for (int phi=0;phi<72;++phi)
00602 
00603                 } // for (int eta=0; eta< (etaBins_-2);++eta)
00604 
00605                 return;
00606         } // FillUnphysicalHEHFBins(std::vector<MonitorElement*> &hh)
00607 
00608 
00609 
00610         // special fill call based on detid -- eventually will need special treatment
00611         void Fill(HcalDetId& id, double val /*=1*/, std::vector<TH2F> &depth)
00612         { 
00613                 // If in HF, need to shift by 1 bin (-1 bin lower in -HF, +1 bin higher in +HF)
00614                 if (id.subdet()==HcalForward)
00615                         depth[id.depth()-1].Fill(id.ieta()<0 ? id.ieta()-1 : id.ieta()+1, id.iphi(), val);
00616                 else 
00617                         depth[id.depth()-1].Fill(id.ieta(),id.iphi(),val);
00618         }
00619 
00620         void Reset(std::vector<TH2F> &depth) 
00621         {
00622                 for (unsigned int d=0;d<depth.size();d++)
00623                         //BUG CAN BE HERE:
00624                         //if(depth[d]) 
00625                         depth[d].Reset();
00626         } // void Reset(void)
00627 
00628         void setup(std::vector<TH2F> &depth, std::string name, std::string units=""){
00629                 std::string unittitle, unitname;
00630                 if (units.empty())
00631                 {
00632                         unitname = units;
00633                         unittitle = "No Units";
00634                 }
00635                 else
00636                 {
00637                         unitname = " " + units;
00638                         unittitle = units;
00639                 }
00640 
00641                 // Push back depth plots
00643                 //1. create first plot  
00644                 depth.push_back(TH2F(("HB HE HF Depth 1 "+name+unitname).c_str(),
00645                         (name+" Depth 1 -- HB HE HF ("+unittitle+")").c_str(),
00646                         85,-42.5,42.5,
00647                         72,0.5,72.5));
00648 
00649                 //2.1 prepare second plot       
00650                 float ybins[73];
00651                 for (int i=0;i<=72;i++) ybins[i]=(float)(i+0.5);
00652                 float xbinsd2[]={-42.5,-41.5,-40.5,-39.5,-38.5,-37.5,-36.5,-35.5,-34.5,-33.5,-32.5,-31.5,-30.5,-29.5,
00653                         -28.5,-27.5,-26.5,-25.5,-24.5,-23.5,-22.5,-21.5,-20.5,-19.5,-18.5,-17.5,-16.5,
00654                         -15.5,-14.5,
00655                         14.5, 15.5,
00656                         16.5,17.5,18.5,19.5,20.5,21.5,22.5,23.5,24.5,25.5,26.5,27.5,28.5,29.5,30.5,
00657                         31.5,32.5,33.5,34.5,35.5,36.5,37.5,38.5,39.5,40.5,41.5,42.5};
00658 
00659                 //2.2 create second plot        
00660                 depth.push_back(TH2F(("HB HE HF Depth 2 "+name+unitname).c_str(),
00661                         (name+" Depth 2 -- HB HE HF ("+unittitle+")").c_str(),
00662                         57, xbinsd2, 72, ybins));
00663 
00664                 //3.1 Set up variable-sized bins for HE depth 3 (MonitorElement also requires phi bins to be entered in array format)
00665                 float xbins[]={-28.5,-27.5,-26.5,-16.5,-15.5,
00666                         15.5,16.5,26.5,27.5,28.5};
00667                 //3.2
00668                 depth.push_back(TH2F(("HE Depth 3 "+name+unitname).c_str(),
00669                         (name+" Depth 3 -- HE ("+unittitle+")").c_str(),
00670                         // Use variable-sized eta bins 
00671                         9, xbins, 72, ybins));
00672 
00673                 //4.1 HO bins are fixed width, but cover a smaller eta range (-15 -> 15)
00674                 depth.push_back(TH2F(("HO Depth 4 "+name+unitname).c_str(),
00675                         (name+" Depth 4 -- HO ("+unittitle+")").c_str(),
00676                         31,-15.5,15.5,
00677                         72,0.5,72.5));
00678 
00679                 for (unsigned int i=0;i<depth.size();++i)
00680                         depth[i].Draw("colz");
00681 
00682                 setBinLabels(depth); // set axis titles, special bins           
00683         }
00684 
00685         void fillOneGain(std::vector<TH2F> &graphData, HcalGains::tAllContWithNames &allContainers, std::string name, int id, std::string units=""){
00686                 setup(graphData, name); 
00687 
00688                 std::stringstream x;
00689                 // Change the titles of each individual histogram
00690                 for (unsigned int d=0;d < graphData.size();++d){
00691                         graphData[d].Reset();
00692                         x << "Gain "<< id  << " for HCAL depth " << d+1;
00693 
00694                         //BUG CAN BE HERE:
00695                         //if (ChannelStatus->depth[d]) 
00696                         graphData[d].SetTitle(x.str().c_str());  // replace "setTitle" with "SetTitle", since you are using TH2F objects instead of MonitorElements
00697                         x.str("");
00698                 }
00699 
00700                 HcalDetId hcal_id;
00701                 int ieta, depth, iphi;
00702 
00703                 //main loop
00704                 // get all containers with names
00705                 //HcalGains::tAllContWithNames allContainers = object().getAllContainers();
00706 
00707                 //ITERATORS AND VALUES:
00708                 HcalGains::tAllContWithNames::const_iterator iter;
00709                 std::vector<HcalGain>::const_iterator contIter;
00710                 float gain = 0.0;
00711 
00712                 //Run trough given id gain:
00713                 int i = id;
00714 
00715                 //run trough all pair containers
00716                 for (iter = allContainers.begin(); iter != allContainers.end(); ++iter){
00717                         //Run trough all values:
00718                         for (contIter = (*iter).second.begin(); contIter != (*iter).second.end(); ++contIter){
00719                                 hcal_id = HcalDetId((uint32_t)(*contIter).rawId());
00720 
00721                                 depth = hcal_id.depth();
00722                                 if (depth<1 || depth>4) 
00723                                         continue;
00724 
00725                                 ieta=hcal_id.ieta();
00726                                 iphi=hcal_id.iphi();
00727 
00728                                 if (hcal_id.subdet() == HcalForward)
00729                                         ieta>0 ? ++ieta : --ieta;
00730 
00731                                 //GET VALUE:
00732                                 gain = (*contIter).getValue(i);
00733                                 //logstatus = log2(1.*channelBits)+1;
00734 
00735                                 //FILLING GOES HERE:
00736                                 graphData[depth-1].Fill(ieta,iphi, gain);
00737 
00738                                 //FOR DEBUGGING:
00739                                 //std::cout << "ieta: " << ieta << "; iphi: " << iphi << "; logstatus: " << logstatus << "; channelBits: " << channelBits<< std::endl;
00740                         }
00741 
00742                 }       
00743         }
00744 
00745 
00746         //FillUnphysicalHEHFBins(graphData);
00747         //return ("kasdasd");
00748 
00749         class ADataRepr //Sample base class for c++ inheritance tutorial
00750         {
00751         public:
00752                 ADataRepr(unsigned int d):m_total(d){};
00753                 unsigned int nr, id;
00754                 std::stringstream filename, rootname, plotname;
00755 
00756                 void fillOneGain(std::vector<TH2F> &graphData, std::string units=""){
00757                         std::stringstream ss("");
00758 
00759                         if (m_total == 1)
00760                                 ss << rootname.str() << " for HCAL depth ";
00761                         else
00762                                 ss << rootname.str() << nr << " for HCAL depth ";
00763 
00764                         setup(graphData, ss.str()); 
00765 
00766                         // Change the titles of each individual histogram
00767                         for (unsigned int d=0;d < graphData.size();++d){
00768                                 graphData[d].Reset();
00769 
00770                                 ss.str("");
00771                                 if (m_total == 1)
00772                                         ss << plotname.str() << " for HCAL depth " << d+1;
00773                                 else
00774                                         ss << plotname.str() << nr << " for HCAL depth " << d+1;
00775 
00776 
00777                                 //BUG CAN BE HERE:
00778                                 //if (ChannelStatus->depth[d]) 
00779                                 graphData[d].SetTitle(ss.str().c_str());  // replace "setTitle" with "SetTitle", since you are using TH2F objects instead of MonitorElements
00780                                 ss.str("");
00781                         }
00782                         //overload this function:
00783                         doFillIn(graphData);
00784 
00785                         FillUnphysicalHEHFBins(graphData);
00786                         
00787                         ss.str("");
00788                         if (m_total == 1)
00789                                 ss << filename.str() << ".png";
00790                         else
00791                                 ss << filename.str() << nr << ".png";
00792                         draw(graphData, ss.str());
00793                         //FOR DEBUGGING:
00794                         //std::cout << "ieta: " << ieta << "; iphi: " << iphi << "; logstatus: " << logstatus << "; channelBits: " << channelBits<< std::endl;
00795                 }
00796 
00797         protected:
00798                 unsigned int m_total;
00799                 HcalDetId hcal_id;
00800                 int ieta, depth, iphi;
00801 
00802 
00803                 virtual void doFillIn(std::vector<TH2F> &graphData) = 0;
00804                 
00805         private:
00806 
00807         void draw(std::vector<TH2F> &graphData, std::string filename) {
00808                 //Drawing...
00809                 // use David's palette
00810                 gStyle->SetPalette(1);
00811                 const Int_t NCont = 999;
00812                 gStyle->SetNumberContours(NCont);
00813                 TCanvas canvas("CC map","CC map",840,369*4);
00814 
00815                 TPad pad1("pad1","pad1", 0.0, 0.75, 1.0, 1.0);
00816                 pad1.Draw();
00817                 TPad pad2("pad2","pad2", 0.0, 0.5, 1.0, 0.75);
00818                 pad2.Draw();
00819                 TPad pad3("pad3","pad3", 0.0, 0.25, 1.0, 0.5);
00820                 pad3.Draw();
00821                 TPad pad4("pad4","pad4", 0.0, 0.0, 1.0, 0.25);
00822                 pad4.Draw();
00823 
00824 
00825                 pad1.cd();
00826                 graphData[0].SetStats(0);
00827                 graphData[0].Draw("colz");
00828 
00829                 pad2.cd();
00830                 graphData[1].SetStats(0);
00831                 graphData[1].Draw("colz");
00832 
00833                 pad3.cd();
00834                 graphData[2].SetStats(0);
00835                 graphData[2].Draw("colz");
00836 
00837                 pad4.cd();
00838                 graphData[3].SetStats(0);
00839                 graphData[3].Draw("colz");
00840 
00841                 canvas.SaveAs(filename.c_str());
00842         }
00843 
00844                 void setup(std::vector<TH2F> &depth, std::string name, std::string units=""){
00845                         std::string unittitle, unitname;
00846                         if (units.empty())
00847                         {
00848                                 unitname = units;
00849                                 unittitle = "No Units";
00850                         }
00851                         else
00852                         {
00853                                 unitname = " " + units;
00854                                 unittitle = units;
00855                         }
00856 
00857                         // Push back depth plots
00859                         //1. create first plot  
00860                         depth.push_back(TH2F(("HB HE HF Depth 1 "+name+unitname).c_str(),
00861                                 (name+" Depth 1 -- HB HE HF ("+unittitle+")").c_str(),
00862                                 85,-42.5,42.5,
00863                                 72,0.5,72.5));
00864 
00865                         //2.1 prepare second plot       
00866                         float ybins[73];
00867                         for (int i=0;i<=72;i++) ybins[i]=(float)(i+0.5);
00868                         float xbinsd2[]={-42.5,-41.5,-40.5,-39.5,-38.5,-37.5,-36.5,-35.5,-34.5,-33.5,-32.5,-31.5,-30.5,-29.5,
00869                                 -28.5,-27.5,-26.5,-25.5,-24.5,-23.5,-22.5,-21.5,-20.5,-19.5,-18.5,-17.5,-16.5,
00870                                 -15.5,-14.5,
00871                                 14.5, 15.5,
00872                                 16.5,17.5,18.5,19.5,20.5,21.5,22.5,23.5,24.5,25.5,26.5,27.5,28.5,29.5,30.5,
00873                                 31.5,32.5,33.5,34.5,35.5,36.5,37.5,38.5,39.5,40.5,41.5,42.5};
00874 
00875                         //2.2 create second plot        
00876                         depth.push_back(TH2F(("HB HE HF Depth 2 "+name+unitname).c_str(),
00877                                 (name+" Depth 2 -- HB HE HF ("+unittitle+")").c_str(),
00878                                 57, xbinsd2, 72, ybins));
00879 
00880                         //3.1 Set up variable-sized bins for HE depth 3 (MonitorElement also requires phi bins to be entered in array format)
00881                         float xbins[]={-28.5,-27.5,-26.5,-16.5,-15.5,
00882                                 15.5,16.5,26.5,27.5,28.5};
00883                         //3.2
00884                         depth.push_back(TH2F(("HE Depth 3 "+name+unitname).c_str(),
00885                                 (name+" Depth 3 -- HE ("+unittitle+")").c_str(),
00886                                 // Use variable-sized eta bins 
00887                                 9, xbins, 72, ybins));
00888 
00889                         //4.1 HO bins are fixed width, but cover a smaller eta range (-15 -> 15)
00890                         depth.push_back(TH2F(("HO Depth 4 "+name+unitname).c_str(),
00891                                 (name+" Depth 4 -- HO ("+unittitle+")").c_str(),
00892                                 31,-15.5,15.5,
00893                                 72,0.5,72.5));
00894 
00895                         for (unsigned int i=0;i<depth.size();++i)
00896                                 depth[i].Draw("colz");
00897 
00898                         setBinLabels(depth); // set axis titles, special bins           
00899                 }
00900 
00901                 //functions for making plot:
00902                 void setBinLabels(std::vector<TH2F> &depth)
00903                 {
00904                         // Set labels for all depth histograms
00905                         for (unsigned int i=0;i<depth.size();++i)
00906                         {
00907                                 depth[i].SetXTitle("i#eta");
00908                                 depth[i].SetYTitle("i#phi");
00909                         }
00910 
00911                         std::stringstream label;
00912 
00913                         // set label on every other bin
00914                         for (int i=-41;i<=-29;i=i+2)
00915                         {
00916                                 label<<i;
00917                                 depth[0].GetXaxis()->SetBinLabel(i+42,label.str().c_str());
00918                                 depth[1].GetXaxis()->SetBinLabel(i+42,label.str().c_str());
00919                                 label.str("");
00920                         }
00921                         depth[0].GetXaxis()->SetBinLabel(14,"-29HE");
00922                         depth[1].GetXaxis()->SetBinLabel(14,"-29HE");
00923 
00924                         // offset by one for HE
00925                         for (int i=-27;i<=27;i=i+2)
00926                         {
00927                                 label<<i;
00928                                 depth[0].GetXaxis()->SetBinLabel(i+43,label.str().c_str());
00929                                 label.str("");
00930                         }
00931                         depth[0].GetXaxis()->SetBinLabel(72,"29HE");
00932                         for (int i=29;i<=41;i=i+2)
00933                         {
00934                                 label<<i;
00935                                 depth[0].GetXaxis()->SetBinLabel(i+44,label.str().c_str());
00936                                 label.str("");
00937                         }
00938                         for (int i=16;i<=28;i=i+2)
00939                         {
00940                                 label<<i-43;
00941                                 depth[1].GetXaxis()->SetBinLabel(i,label.str().c_str());
00942                                 label.str("");
00943                         }
00944                         depth[1].GetXaxis()->SetBinLabel(29,"NULL");
00945                         for (int i=15;i<=27;i=i+2)
00946                         {
00947                                 label<<i;
00948                                 depth[1].GetXaxis()->SetBinLabel(i+15,label.str().c_str());
00949                                 label.str("");
00950                         }
00951 
00952                         depth[1].GetXaxis()->SetBinLabel(44,"29HE");
00953                         for (int i=29;i<=41;i=i+2)
00954                         {
00955                                 label<<i;
00956                                 depth[1].GetXaxis()->SetBinLabel(i+16,label.str().c_str());
00957                                 label.str("");
00958                         }
00959 
00960                         // HE depth 3 labels;
00961                         depth[2].GetXaxis()->SetBinLabel(1,"-28");
00962                         depth[2].GetXaxis()->SetBinLabel(2,"-27");
00963                         depth[2].GetXaxis()->SetBinLabel(3,"Null");
00964                         depth[2].GetXaxis()->SetBinLabel(4,"-16");
00965                         depth[2].GetXaxis()->SetBinLabel(5,"Null");
00966                         depth[2].GetXaxis()->SetBinLabel(6,"16");
00967                         depth[2].GetXaxis()->SetBinLabel(7,"Null");
00968                         depth[2].GetXaxis()->SetBinLabel(8,"27");
00969                         depth[2].GetXaxis()->SetBinLabel(9,"28");
00970                 }
00971 
00972         };
00973 }
00974 #endif