CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/DQMOffline/CalibCalo/src/DQMHcalPhiSymAlCaReco.cc

Go to the documentation of this file.
00001 /*
00002  * \file DQMHcalPhiSymAlCaReco.cc
00003  *
00004  * \author Olga Kodolova
00005  * 
00006  * $Date: 2010/10/15 22:44:31 $
00007  * $Revision: 1.17 $
00008  *
00009  *
00010  * Description: Monitoring of Phi Symmetry Calibration Stream  
00011 */
00012 
00013 #include "FWCore/Framework/interface/Event.h"
00014 #include "FWCore/Framework/interface/EventSetup.h"
00015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00016 #include "FWCore/ServiceRegistry/interface/Service.h"
00017 
00018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00019 
00020 // DQM include files
00021 
00022 #include "DQMServices/Core/interface/MonitorElement.h"
00023 
00024 // work on collections
00025 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
00026 
00027 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00028 #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
00029 #include "DataFormats/HcalDetId/interface/HcalSubdetector.h"
00030 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
00031 
00032 #include "DataFormats/FEDRawData/interface/FEDRawDataCollection.h"
00033 #include "DataFormats/FEDRawData/interface/FEDRawData.h"
00034 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
00035 #include "DataFormats/FEDRawData/interface/FEDHeader.h"
00036 #include "EventFilter/HcalRawToDigi/interface/HcalHTRData.h"
00037 #include "EventFilter/HcalRawToDigi/interface/HcalDCCHeader.h"
00038 
00039 #include "DQMServices/Core/interface/DQMStore.h"
00040 #include "DQMOffline/CalibCalo/src/DQMHcalPhiSymAlCaReco.h"
00041 
00042 using namespace std;
00043 using namespace edm;
00044 
00045 // ******************************************
00046 // constructors
00047 // *****************************************
00048 
00049 DQMHcalPhiSymAlCaReco::DQMHcalPhiSymAlCaReco( const edm::ParameterSet& ps ) :
00050 eventCounter_(0)
00051 {
00052   dbe_ = Service<DQMStore>().operator->();
00053   //
00054   // Input from configurator file 
00055   //
00056   folderName_   = ps.getUntrackedParameter<string>("FolderName","ALCAStreamHcalPhiSym");
00057   
00058   hbherecoMB    = ps.getParameter<edm::InputTag>("hbheInputMB");
00059   horecoMB      = ps.getParameter<edm::InputTag>("hoInputMB");
00060   hfrecoMB      = ps.getParameter<edm::InputTag>("hfInputMB");
00061   
00062   hbherecoNoise = ps.getParameter<edm::InputTag>("hbheInputNoise");
00063   horecoNoise   = ps.getParameter<edm::InputTag>("hoInputNoise");
00064   hfrecoNoise   = ps.getParameter<edm::InputTag>("hfInputNoise");
00065 
00066   rawInLabel_   = ps.getParameter<edm::InputTag>("rawInputLabel");
00067 
00068   period_       = ps.getParameter<unsigned int>("period") ;
00069 
00070   saveToFile_   = ps.getUntrackedParameter<bool>("SaveToFile",false);
00071   fileName_     = ps.getUntrackedParameter<string>("FileName","MonitorAlCaHcalPhiSym.root");
00072 
00073   // histogram parameters
00074 
00075   // Distribution of rechits in iPhi, iEta 
00076   hiDistr_y_nbin_ = ps.getUntrackedParameter<int>("hiDistr_y_nbin",72);
00077   hiDistr_y_min_  =  ps.getUntrackedParameter<double>("hiDistr_y_min",0.5);
00078   hiDistr_y_max_  =  ps.getUntrackedParameter<double>("hiDistr_y_max",72.5);
00079   hiDistr_x_nbin_ = ps.getUntrackedParameter<int>("hiDistr_x_nbin",41);
00080   hiDistr_x_min_  =  ps.getUntrackedParameter<double>("hiDistr_x_min",0.5);
00081   hiDistr_x_max_  =  ps.getUntrackedParameter<double>("hiDistr_x_max",41.5);
00082   // Check for NZS
00083   hiDistr_r_nbin_ = ps.getUntrackedParameter<int>("hiDistr_r_nbin",100);
00084   ihbhe_size_     = ps.getUntrackedParameter<double>("ihbhe_size_",5184.);
00085   ihf_size_       = ps.getUntrackedParameter<double>("ihf_size_",1728.);
00086     
00087 }
00088 
00089 DQMHcalPhiSymAlCaReco::~DQMHcalPhiSymAlCaReco()
00090 {}
00091 
00092 //--------------------------------------------------------
00093 void DQMHcalPhiSymAlCaReco::beginJob(){
00094  
00095   // create and cd into new folder
00096   dbe_->setCurrentFolder(folderName_);
00097 
00098   eventCounter_ = 0;
00099 
00100   hFEDsize = dbe_->book1D("hFEDsize","HCAL FED size (kB)",200,-0.5,20.5);
00101   hFEDsize->setAxisTitle("kB",1);
00102 
00103   hHcalIsZS = dbe_->book1D("hHcalIsZS", "Hcal Is ZS", 4, -1.5, 2.5);
00104   hHcalIsZS->setBinLabel(2, "NZS");
00105   hHcalIsZS->setBinLabel(3, "ZS");
00106 
00107   char hname[50];
00108   sprintf(hname, "L1 Event Number %% %i", period_);
00109   hL1Id = dbe_->book1D("hL1Id", hname,4200,-99.5,4099.5);
00110   hL1Id->setAxisTitle(hname);
00111   
00112 
00113   // book some histograms 1D
00114   double xmin = 0.1;
00115   double xmax = 1.1;
00116   hiDistrHBHEsize1D_ =
00117     dbe_->book1D("DistrHBHEsize","Size of HBHE Collection",
00118                   hiDistr_r_nbin_,
00119                   xmin,
00120                   xmax
00121                 );
00122   hiDistrHFsize1D_ =
00123     dbe_->book1D("DistrHFsize","Size of HF Collection",
00124                   hiDistr_r_nbin_,
00125                   xmin,
00126                   xmax
00127                 );
00128 
00129 
00130   // First moment
00131   hiDistrMBPl2D_ = 
00132     dbe_->book2D("MBdepthPl1", "iphi- +ieta signal distribution at depth1",
00133                  hiDistr_x_nbin_, 
00134                  hiDistr_x_min_,
00135                  hiDistr_x_max_,
00136                  hiDistr_y_nbin_, 
00137                  hiDistr_y_min_,
00138                  hiDistr_y_max_
00139                  );
00140 
00141   hiDistrMBPl2D_->setAxisTitle("i#phi ", 2);
00142   hiDistrMBPl2D_->setAxisTitle("i#eta ", 1);
00143 
00144 
00145   hiDistrNoisePl2D_ = 
00146     dbe_->book2D("NoisedepthPl1", "iphi-ieta noise distribution at depth1",
00147                  hiDistr_x_nbin_+1, 
00148                  hiDistr_x_min_-1.,
00149                  hiDistr_x_max_,
00150                  hiDistr_y_nbin_+1, 
00151                  hiDistr_y_min_-1.,
00152                  hiDistr_y_max_
00153                  );
00154 
00155   hiDistrNoisePl2D_->setAxisTitle("i#phi ", 2);
00156   hiDistrNoisePl2D_->setAxisTitle("i#eta ", 1);
00157 // Second moment
00158   hiDistrMB2Pl2D_ =
00159     dbe_->book2D("MB2depthPl1", "iphi- +ieta signal distribution at depth1",
00160                  hiDistr_x_nbin_,
00161                  hiDistr_x_min_,
00162                  hiDistr_x_max_,
00163                  hiDistr_y_nbin_,
00164                  hiDistr_y_min_,
00165                  hiDistr_y_max_
00166                  );
00167 
00168   hiDistrMB2Pl2D_->setAxisTitle("i#phi ", 2);
00169   hiDistrMB2Pl2D_->setAxisTitle("i#eta ", 1);
00170 
00171 
00172   hiDistrNoise2Pl2D_ =
00173     dbe_->book2D("Noise2depthPl1", "iphi-ieta noise distribution at depth1",
00174                  hiDistr_x_nbin_,
00175                  hiDistr_x_min_,
00176                  hiDistr_x_max_,
00177                  hiDistr_y_nbin_,
00178                  hiDistr_y_min_,
00179                  hiDistr_y_max_
00180                  );
00181 
00182   hiDistrNoise2Pl2D_->setAxisTitle("i#phi ", 2);
00183   hiDistrNoise2Pl2D_->setAxisTitle("i#eta ", 1);
00184 
00185 // Variance
00186   hiDistrVarMBPl2D_ =
00187     dbe_->book2D("VarMBdepthPl1", "iphi- +ieta signal distribution at depth1",
00188                  hiDistr_x_nbin_,
00189                  hiDistr_x_min_,
00190                  hiDistr_x_max_,
00191                  hiDistr_y_nbin_,
00192                  hiDistr_y_min_,
00193                  hiDistr_y_max_
00194                  );
00195 
00196   hiDistrVarMBPl2D_->setAxisTitle("i#phi ", 2);
00197   hiDistrVarMBPl2D_->setAxisTitle("i#eta ", 1);
00198 
00199 
00200   hiDistrVarNoisePl2D_ =
00201     dbe_->book2D("VarNoisedepthPl1", "iphi-ieta noise distribution at depth1",
00202                  hiDistr_x_nbin_,
00203                  hiDistr_x_min_,
00204                  hiDistr_x_max_,
00205                  hiDistr_y_nbin_,
00206                  hiDistr_y_min_,
00207                  hiDistr_y_max_
00208                  );
00209 
00210   hiDistrVarNoisePl2D_->setAxisTitle("i#phi ", 2);
00211   hiDistrVarNoisePl2D_->setAxisTitle("i#eta ", 1);
00212 
00213 //==================================================================================
00214 // First moment
00215   hiDistrMBMin2D_ = 
00216     dbe_->book2D("MBdepthMin1", "iphi- +ieta signal distribution at depth1",
00217                  hiDistr_x_nbin_, 
00218                  hiDistr_x_min_,
00219                  hiDistr_x_max_,
00220                  hiDistr_y_nbin_, 
00221                  hiDistr_y_min_,
00222                  hiDistr_y_max_
00223                  );
00224 
00225   hiDistrMBMin2D_->setAxisTitle("i#phi ", 2);
00226   hiDistrMBMin2D_->setAxisTitle("i#eta ", 1);
00227 
00228 
00229   hiDistrNoiseMin2D_ = 
00230     dbe_->book2D("NoisedepthMin1", "iphi-ieta noise distribution at depth1",
00231                  hiDistr_x_nbin_, 
00232                  hiDistr_x_min_,
00233                  hiDistr_x_max_,
00234                  hiDistr_y_nbin_, 
00235                  hiDistr_y_min_,
00236                  hiDistr_y_max_
00237                  );
00238 
00239   hiDistrNoiseMin2D_->setAxisTitle("i#phi ", 2);
00240   hiDistrNoiseMin2D_->setAxisTitle("i#eta ", 1);
00241 // Second moment
00242   hiDistrMB2Min2D_ =
00243     dbe_->book2D("MB2depthMin1", "iphi- +ieta signal distribution at depth1",
00244                  hiDistr_x_nbin_,
00245                  hiDistr_x_min_,
00246                  hiDistr_x_max_,
00247                  hiDistr_y_nbin_,
00248                  hiDistr_y_min_,
00249                  hiDistr_y_max_
00250                  );
00251 
00252   hiDistrMB2Min2D_->setAxisTitle("i#phi ", 2);
00253   hiDistrMB2Min2D_->setAxisTitle("i#eta ", 1);
00254 
00255 
00256   hiDistrNoise2Min2D_ =
00257     dbe_->book2D("Noise2depthMin1", "iphi-ieta noise distribution at depth1",
00258                  hiDistr_x_nbin_,
00259                  hiDistr_x_min_,
00260                  hiDistr_x_max_,
00261                  hiDistr_y_nbin_,
00262                  hiDistr_y_min_,
00263                  hiDistr_y_max_
00264                  );
00265 
00266   hiDistrNoise2Min2D_->setAxisTitle("i#phi ", 2);
00267   hiDistrNoise2Min2D_->setAxisTitle("i#eta ", 1);
00268 
00269 // Variance
00270   hiDistrVarMBMin2D_ =
00271     dbe_->book2D("VarMBdepthMin1", "iphi- +ieta signal distribution at depth1",
00272                  hiDistr_x_nbin_,
00273                  hiDistr_x_min_,
00274                  hiDistr_x_max_,
00275                  hiDistr_y_nbin_,
00276                  hiDistr_y_min_,
00277                  hiDistr_y_max_
00278                  );
00279 
00280   hiDistrVarMBMin2D_->setAxisTitle("i#phi ", 2);
00281   hiDistrVarMBMin2D_->setAxisTitle("i#eta ", 1);
00282 
00283 
00284   hiDistrVarNoiseMin2D_ =
00285     dbe_->book2D("VarNoisedepthMin1", "iphi-ieta noise distribution at depth1",
00286                  hiDistr_x_nbin_,
00287                  hiDistr_x_min_,
00288                  hiDistr_x_max_,
00289                  hiDistr_y_nbin_,
00290                  hiDistr_y_min_,
00291                  hiDistr_y_max_
00292                  );
00293 
00294   hiDistrVarNoiseMin2D_->setAxisTitle("i#phi ", 2);
00295   hiDistrVarNoiseMin2D_->setAxisTitle("i#eta ", 1);
00296 
00297 
00298 }
00299 
00300 //--------------------------------------------------------
00301 void DQMHcalPhiSymAlCaReco::beginRun(const edm::Run& r, const EventSetup& context) {
00302 //   eventCounter_ = 0;
00303 }
00304 
00305 //--------------------------------------------------------
00306 void DQMHcalPhiSymAlCaReco::beginLuminosityBlock(const LuminosityBlock& lumiSeg, 
00307      const EventSetup& context) {
00308   
00309 }
00310 
00311 //-------------------------------------------------------------
00312 
00313 void DQMHcalPhiSymAlCaReco::analyze(const Event& iEvent, 
00314                                const EventSetup& iSetup ){  
00315  
00316 
00317   eventCounter_++;
00318  
00319   edm::Handle<FEDRawDataCollection> rawIn;
00320   iEvent.getByLabel(rawInLabel_,rawIn);
00321 
00322   if(!rawIn.isValid()){
00323      LogDebug("") << "HcalCalibAlgos: Error! can't get hbhe product!" << std::endl;
00324       return ;
00325   }
00326 
00327   //get HCAL FEDs:
00328   std::vector<int> selFEDs;
00329   for (int i=FEDNumbering::MINHCALFEDID; i<=FEDNumbering::MAXHCALFEDID; i++)
00330     {
00331       selFEDs.push_back(i);
00332     }
00333 
00334   //  std::cout<<" Size of FED "<<selFEDs.size()<<std::endl;
00335 
00336   const FEDRawDataCollection *rdc=rawIn.product(); 
00337 
00338   bool hcalIsZS = false ;
00339   for (unsigned int k=0; k<selFEDs.size(); k++)
00340     {
00341       const FEDRawData & fedData = rdc->FEDData(selFEDs[k]);
00342       // std::cout<<fedData.size()*std::pow(1024.,-1)<<std::endl;
00343       hFEDsize->Fill(fedData.size()*std::pow(1024.,-1),1);
00344       
00345       // get HCAL DCC Header for each FEDRawData
00346       const HcalDCCHeader* dccHeader=(const HcalDCCHeader*)(fedData.data());
00347 
00348       if (!dccHeader) continue;
00349 
00350       // walk through the HTR data...
00351       HcalHTRData htr;
00352 
00353       int nspigot =0; 
00354       for (int spigot=0; spigot<HcalDCCHeader::SPIGOT_COUNT; spigot++) {    
00355        nspigot++;
00356 
00357         if (!dccHeader->getSpigotPresent(spigot)) continue;
00358 
00359         // Load the given decoder with the pointer and length from this spigot.
00360         dccHeader->getSpigotData(spigot,htr, fedData.size()); 
00361 
00362         if(k != 20 && nspigot !=14 ) {      
00363         if ( !htr.isUnsuppressed() ) { hcalIsZS = true; }
00364         }
00365       } 
00366     
00367     } // loop over HcalFEDs
00368   
00369   hHcalIsZS->Fill( hcalIsZS );
00370 
00371   // get Trigger FED-Id
00372   const FEDRawData& fedData = rdc->FEDData(FEDNumbering::MINTriggerGTPFEDID) ;
00373   FEDHeader header(fedData.data()) ;
00374 
00376   hL1Id->Fill( (header.lvl1ID())%period_ );
00377  
00378   edm::Handle<HBHERecHitCollection> hbheNS;
00379   iEvent.getByLabel(hbherecoNoise, hbheNS);
00380   
00381   if(!hbheNS.isValid()){
00382     LogDebug("") << "HcalCalibAlgos: Error! can't get hbhe product!" << std::endl;
00383     return ;
00384   }
00385   
00386   edm::Handle<HBHERecHitCollection> hbheMB;
00387   iEvent.getByLabel(hbherecoMB, hbheMB);
00388   
00389   if(!hbheMB.isValid()){
00390     LogDebug("") << "HcalCalibAlgos: Error! can't get hbhe product!" << std::endl;
00391     return ;
00392   }
00393   
00394   edm::Handle<HFRecHitCollection> hfNS;
00395   iEvent.getByLabel(hfrecoNoise, hfNS);
00396   
00397   if(!hfNS.isValid()){
00398     LogDebug("") << "HcalCalibAlgos: Error! can't get hbhe product!" << std::endl;
00399     return ;
00400   }
00401   
00402   edm::Handle<HFRecHitCollection> hfMB;
00403   iEvent.getByLabel(hfrecoMB, hfMB);
00404   
00405   if(!hfMB.isValid()){
00406     LogDebug("") << "HcalCalibAlgos: Error! can't get hbhe product!" << std::endl;
00407     return ;
00408   }
00409   
00410   const HBHERecHitCollection HithbheNS = *(hbheNS.product());
00411   
00412   hiDistrHBHEsize1D_->Fill(HithbheNS.size()/ihbhe_size_);
00413   
00414      
00415   for(HBHERecHitCollection::const_iterator hbheItr=HithbheNS.begin(); hbheItr!=HithbheNS.end(); hbheItr++)
00416         {
00417                  DetId id = (*hbheItr).detid(); 
00418                  HcalDetId hid=HcalDetId(id);
00419                  
00420                  if(hid.depth() == 1) {
00421                  if( hid.ieta() > 0 ) {
00422                  hiDistrNoisePl2D_->Fill(hid.ieta(),hid.iphi(),hbheItr->energy());
00423                  hiDistrNoise2Pl2D_->Fill(hid.ieta(),hid.iphi(),hbheItr->energy()*hbheItr->energy());
00424                  } else {
00425                  hiDistrNoiseMin2D_->Fill(fabs(hid.ieta()),hid.iphi(),hbheItr->energy());
00426                  hiDistrNoise2Min2D_->Fill(fabs(hid.ieta()),hid.iphi(),hbheItr->energy()*hbheItr->energy());
00427                  }
00428                  }
00429         }
00430 
00431   const HBHERecHitCollection HithbheMB = *(hbheMB.product());
00432 
00433   for(HBHERecHitCollection::const_iterator hbheItr=HithbheMB.begin(); hbheItr!=HithbheMB.end(); hbheItr++)
00434         {
00435                  DetId id = (*hbheItr).detid(); 
00436                  HcalDetId hid=HcalDetId(id);
00437                  
00438                  if(hid.depth() == 1) {
00439                  if( hid.ieta() > 0 ) {
00440                  hiDistrMBPl2D_->Fill(hid.ieta(),hid.iphi(),hbheItr->energy());
00441                  hiDistrMB2Pl2D_->Fill(hid.ieta(),hid.iphi(),hbheItr->energy()*hbheItr->energy());
00442                  } else {
00443                  hiDistrMBMin2D_->Fill(fabs(hid.ieta()),hid.iphi(),hbheItr->energy());
00444                  hiDistrMB2Min2D_->Fill(fabs(hid.ieta()),hid.iphi(),hbheItr->energy()*hbheItr->energy());
00445                  }
00446                  }
00447 
00448         }
00449 
00450   const HFRecHitCollection HithfNS = *(hfNS.product());
00451 
00452   hiDistrHFsize1D_->Fill(HithfNS.size()/ihf_size_);
00453   
00454   for(HFRecHitCollection::const_iterator hbheItr=HithfNS.begin(); hbheItr!=HithfNS.end(); hbheItr++)
00455         {
00456         
00457                  DetId id = (*hbheItr).detid(); 
00458                  HcalDetId hid=HcalDetId(id);
00459                  
00460                  if(hid.depth() == 1) {
00461                  if( hid.ieta() > 0 ) {
00462                  hiDistrNoisePl2D_->Fill(hid.ieta(),hid.iphi(),hbheItr->energy());
00463                  hiDistrNoise2Pl2D_->Fill(hid.ieta(),hid.iphi(),hbheItr->energy()*hbheItr->energy());
00464                  } else {
00465                  hiDistrNoiseMin2D_->Fill(fabs(hid.ieta()),hid.iphi(),hbheItr->energy());
00466                  hiDistrNoise2Min2D_->Fill(fabs(hid.ieta()),hid.iphi(),hbheItr->energy()*hbheItr->energy());
00467                  }
00468                  }
00469         
00470         }
00471 
00472   const HFRecHitCollection HithfMB = *(hfMB.product());
00473  
00474   for(HFRecHitCollection::const_iterator hbheItr=HithfMB.begin(); hbheItr!=HithfMB.end(); hbheItr++)
00475         {
00476                  DetId id = (*hbheItr).detid(); 
00477                  HcalDetId hid=HcalDetId(id);
00478                  
00479                  if(hid.depth() == 1) {
00480                  if( hid.ieta() > 0 ) {
00481                  hiDistrMBPl2D_->Fill(hid.ieta(),hid.iphi(),hbheItr->energy());
00482                  hiDistrMB2Pl2D_->Fill(hid.ieta(),hid.iphi(),hbheItr->energy()*hbheItr->energy());
00483                  } else {
00484                  hiDistrMBMin2D_->Fill(fabs(hid.ieta()),hid.iphi(),hbheItr->energy());
00485                  hiDistrMB2Min2D_->Fill(fabs(hid.ieta()),hid.iphi(),hbheItr->energy()*hbheItr->energy());
00486                  }
00487                  }
00488         }       
00489         
00490         
00491 } //analyze
00492 
00493 
00494 
00495 
00496 //--------------------------------------------------------
00497 void DQMHcalPhiSymAlCaReco::endLuminosityBlock(const LuminosityBlock& lumiSeg, 
00498                                           const EventSetup& context) {
00499 }
00500 //--------------------------------------------------------
00501 void DQMHcalPhiSymAlCaReco::endRun(const Run& r, const EventSetup& context){
00502 // Keep Variances
00503   if(eventCounter_ > 0) {
00504   for(int k=0; k<=hiDistr_x_nbin_;k++)
00505   {
00506     for(int j=0; j<=hiDistr_y_nbin_;j++)
00507     {
00508 // First moment
00509        float cc1=hiDistrMBPl2D_->getBinContent(k,j);
00510        cc1 = cc1 * 1./eventCounter_;
00511        float cc2=hiDistrNoisePl2D_->getBinContent(k,j);
00512        cc2 = cc2 * 1./eventCounter_;
00513        float cc3=hiDistrMBMin2D_->getBinContent(k,j);
00514        cc3 = cc3 * 1./eventCounter_;
00515        float cc4=hiDistrNoiseMin2D_->getBinContent(k,j);
00516        cc4 = cc4 * 1./eventCounter_;
00517 // Second moment
00518        float cc11=hiDistrMB2Pl2D_->getBinContent(k,j);
00519        cc11 = cc11 * 1./eventCounter_;
00520        hiDistrVarMBPl2D_->setBinContent(k,j,cc11-cc1*cc1);
00521        float cc22=hiDistrNoise2Pl2D_->getBinContent(k,j);
00522        cc22 = cc22 * 1./eventCounter_;
00523        hiDistrVarNoisePl2D_->setBinContent(k,j,cc22-cc2*cc2);
00524        float cc33=hiDistrMB2Min2D_->getBinContent(k,j);
00525        cc33 = cc33 * 1./eventCounter_;
00526        hiDistrVarMBMin2D_->setBinContent(k,j,cc33-cc3*cc3);
00527        float cc44=hiDistrNoise2Min2D_->getBinContent(k,j);
00528        cc44 = cc44 * 1./eventCounter_;
00529        hiDistrVarNoiseMin2D_->setBinContent(k,j,cc44-cc4*cc4);
00530     }     
00531   }
00532   }
00533 }
00534 //--------------------------------------------------------
00535 void DQMHcalPhiSymAlCaReco::endJob(){
00536   if (saveToFile_) {
00537 
00538   for(int k=0; k<=hiDistr_x_nbin_;k++)
00539   {
00540     for(int j=0; j<=hiDistr_y_nbin_;j++)
00541     {
00542 // First moment
00543        float cc1=hiDistrMBPl2D_->getBinContent(k,j);
00544        cc1 = cc1 * 1./eventCounter_;
00545        hiDistrMBPl2D_->setBinContent(k,j,cc1); 
00546        float cc2=hiDistrNoisePl2D_->getBinContent(k,j);
00547        cc2 = cc2 * 1./eventCounter_;
00548        hiDistrNoisePl2D_->setBinContent(k,j,cc2);
00549        float cc3=hiDistrMBMin2D_->getBinContent(k,j);
00550        cc3 = cc3 * 1./eventCounter_;
00551        hiDistrMBMin2D_->setBinContent(k,j,cc3);
00552        float cc4=hiDistrNoiseMin2D_->getBinContent(k,j);
00553        cc4 = cc4 * 1./eventCounter_;
00554        hiDistrNoiseMin2D_->setBinContent(k,j,cc4);
00555 // Second moment
00556        float cc11=hiDistrMB2Pl2D_->getBinContent(k,j);
00557        cc11 = cc11 * 1./eventCounter_;
00558        hiDistrMB2Pl2D_->setBinContent(k,j,cc11);
00559        hiDistrVarMBPl2D_->setBinContent(k,j,cc11-cc1*cc1);
00560        float cc22=hiDistrNoise2Pl2D_->getBinContent(k,j);
00561        cc22 = cc22 * 1./eventCounter_;
00562        hiDistrNoise2Pl2D_->setBinContent(k,j,cc22);
00563        hiDistrVarNoisePl2D_->setBinContent(k,j,cc22-cc2*cc2);
00564        float cc33=hiDistrMB2Min2D_->getBinContent(k,j);
00565        cc33 = cc33 * 1./eventCounter_;
00566        hiDistrMB2Min2D_->setBinContent(k,j,cc33);
00567        hiDistrVarMBMin2D_->setBinContent(k,j,cc33-cc3*cc3);
00568        float cc44=hiDistrNoise2Min2D_->getBinContent(k,j);
00569        cc44 = cc44 * 1./eventCounter_;
00570        hiDistrNoise2Min2D_->setBinContent(k,j,cc44);
00571        hiDistrVarNoiseMin2D_->setBinContent(k,j,cc44-cc4*cc4);
00572     }
00573   }
00574      dbe_->save(fileName_);
00575   }
00576 }
00577 
00578