CMS 3D CMS Logo

HcalPedestalMonitor.cc

Go to the documentation of this file.
00001 #include "DQM/HcalMonitorTasks/interface/HcalPedestalMonitor.h"
00002 #include "DQMServices/Core/interface/DQMStore.h"
00003 #include "DQMServices/Core/interface/MonitorElement.h"
00004 
00005 HcalPedestalMonitor::HcalPedestalMonitor() 
00006 { 
00007   shape_=NULL; 
00008 } // constructor
00009 
00010 
00011 HcalPedestalMonitor::~HcalPedestalMonitor() 
00012 {
00013   // Do we need to delete all pointers here?  If not, will he have a memory leak?  Does this even get explicitly called?  cout statements placed here didn't seem to work
00014   
00015 } // destructor
00016 
00017 
00018 void HcalPedestalMonitor::reset(){}
00019 
00020 
00021 void HcalPedestalMonitor::clearME()
00022 {
00023   // remove monitor elements.  Is this necessary?
00024   if(m_dbe)
00025     {
00026       m_dbe->setCurrentFolder(baseFolder_);
00027       m_dbe->removeContents();
00028     }
00029   return;
00030 } // void HcalPedestalMonitor::clearME();
00031 
00032 
00033 void HcalPedestalMonitor::setup(const edm::ParameterSet& ps, DQMStore* dbe)
00034 {
00035   HcalBaseMonitor::setup(ps,dbe);
00036 
00037   if (showTiming)
00038     {
00039       cpu_timer.reset(); cpu_timer.start();
00040     }
00041 
00042   if (fVerbosity)
00043     cout <<"<HcalPedestalMonitor::setup>  Setting up histograms"<<endl;
00044 
00045   baseFolder_ = rootFolder_+"PedestalMonitor_Hcal";
00046 
00047   // Pedestal Monitor - specific cfg variables
00048 
00049   // set number of ped events needed for pedestal computation to be performed
00050 
00051   doFCpeds_ = ps.getUntrackedParameter<bool>("PedestalMonitor_pedestalsInFC", true);
00052 
00053   minEntriesPerPed_ = ps.getUntrackedParameter<unsigned int>("PedestalMonitor_minEntriesPerPed",1);
00054 
00055   // set expected pedestal mean, width (in ADC)
00056   nominalPedMeanInADC_ = ps.getUntrackedParameter<double>("PedestalMonitor_nominalPedMeanInADC",3);
00057   nominalPedWidthInADC_ = ps.getUntrackedParameter<double>("PedestalMonitor_nominalPedWidthInADC",1);
00058 
00059   // Set error limits that will cause problem histograms to be filled
00060   maxPedMeanDiffADC_ = ps.getUntrackedParameter<double>("PedesstalMonitor_maxPedMeanDiffADC",1.);
00061   maxPedWidthDiffADC_ = ps.getUntrackedParameter<double>("PedestalMonitor_maxPedWidthDiffADC",1.);
00062 
00063   pedmon_minErrorFlag_ = ps.getUntrackedParameter<double>("PedestalMonitor_minErrorFlag", minErrorFlag_);
00064   pedmon_checkNevents_ = ps.getUntrackedParameter<int>("PedestalMonitor_checkNevents", checkNevents_);
00065   // set bins over which pedestals will be computed
00066   startingTimeSlice_ = ps.getUntrackedParameter<int>("PedestalMonitor_startingTimeSlice",0);
00067   endingTimeSlice_   = ps.getUntrackedParameter<int>("PedestalMonitor_endingTimeSlice"  ,1); 
00068 
00069   ievt_=0;
00070 
00071   if ( m_dbe ) 
00072     {
00073       m_dbe->setCurrentFolder(baseFolder_);
00074       meEVT_ = m_dbe->bookInt("Pedestal Task Event Number");
00075       meEVT_->Fill(ievt_);
00076 
00077       // MeanMap, RMSMap values are in ADC, because they show the cells that can create problem cell errors 
00078       setupDepthHists2D(MeanMapByDepth,"Pedestal Mean Map", "ADC");
00079       setupDepthHists2D(RMSMapByDepth, "Pedestal RMS Map", "ADC");
00080       
00081       ProblemPedestals=m_dbe->book2D(" ProblemPedestals",
00082                                      " Problem Pedestal Rate for all HCAL",
00083                                      etaBins_,etaMin_,etaMax_,
00084                                      phiBins_,phiMin_,phiMax_);
00085       ProblemPedestals->setAxisTitle("i#eta",1);
00086       ProblemPedestals->setAxisTitle("i#phi",2);
00087       
00088       // Overall Problem plot appears in main directory; plots by depth appear in subdirectory
00089       m_dbe->setCurrentFolder(baseFolder_+"/problem_pedestals");
00090 
00091       setupDepthHists2D(ProblemPedestalsByDepth, " Problem Pedestal Rate","");
00092 
00093       m_dbe->setCurrentFolder(baseFolder_+"/adc/raw");
00094       setupDepthHists2D(rawADCPedestalMean, "Pedestal Values Map","ADC");
00095       setupDepthHists2D( rawADCPedestalRMS, "Pedestal Widths Map","ADC");
00096       setupDepthHists1D(rawADCPedestalMean_1D, "1D Pedestal Values",
00097                         "ADC",0,10,200);
00098       setupDepthHists1D(rawADCPedestalRMS_1D, "1D Pedestal Widths",
00099                         "ADC",0,10,200);
00100       m_dbe->setCurrentFolder(baseFolder_+"/adc/subtracted__beta_testing");
00101       setupDepthHists2D(subADCPedestalMean, "Subtracted Pedestal Values Map",
00102                         "ADC");
00103       setupDepthHists2D(subADCPedestalRMS, "Subtracted Pedestal Widths Map",
00104                         "ADC");
00105       setupDepthHists1D(subADCPedestalMean_1D, "1D Subtracted Pedestal Values",
00106                         "ADC",-5,5,200);
00107       setupDepthHists1D(subADCPedestalRMS_1D, "1D Subtracted Pedestal Widths",
00108                         "ADC",-5,5,200);
00109 
00110       m_dbe->setCurrentFolder(baseFolder_+"/fc/raw");
00111       setupDepthHists2D(rawFCPedestalMean, "Pedestal Values Map",
00112                         "fC");
00113       setupDepthHists2D(rawFCPedestalRMS, "Pedestal Widths Map",
00114                         "fC");
00115       setupDepthHists1D(rawFCPedestalMean_1D, "1D Pedestal Values",
00116                         "fC",-5,15,200);
00117       setupDepthHists1D(rawFCPedestalRMS_1D, "1D Pedestal Widths",
00118                         "fC",0,10,200);
00119       m_dbe->setCurrentFolder(baseFolder_+"/fc/subtracted__beta_testing");
00120       setupDepthHists2D(subFCPedestalMean, "Subtracted Pedestal Values Map",
00121                         "fC");
00122       setupDepthHists2D(subFCPedestalRMS, "Subtracted Pedestal Widths Map",
00123                         "fC");
00124       setupDepthHists1D(subFCPedestalMean_1D, "1D Subtracted Pedestal Values",
00125                         "fC",-10,10,200);
00126       setupDepthHists1D(subFCPedestalRMS_1D, "1D Subtracted Pedestal Widths",
00127                         "fC",-5,5,200);
00128 
00129       m_dbe->setCurrentFolder(baseFolder_+"/reference_pedestals/adc");
00130       setupDepthHists2D(ADC_PedestalFromDB, ADC_PedestalFromDBByDepth, 
00131                         "Pedestal Values from DataBase","ADC");
00132       setupDepthHists2D(ADC_WidthFromDB, ADC_WidthFromDBByDepth, 
00133                         "Pedestal Widths from DataBase","ADC");
00134 
00135       m_dbe->setCurrentFolder(baseFolder_+"/reference_pedestals/fc");
00136       setupDepthHists2D(fC_PedestalFromDB, fC_PedestalFromDBByDepth, 
00137                         "Pedestal Values from DataBase","fC");
00138                       
00139       setupDepthHists2D(fC_WidthFromDB, fC_WidthFromDBByDepth, 
00140                         "Pedestal Widths from DataBase","fC");
00141 
00142       // initialize all counters to 0
00143       for (int eta=0;eta<ETABINS;++eta)
00144         {
00145           for (int phi=0;phi<PHIBINS;++phi)
00146             {
00147               for (int depth=0;depth<6;++depth)
00148                 {
00149                   pedcounts[eta][phi][depth]=0;
00150                   rawpedsum[eta][phi][depth]=0;
00151                   rawpedsum2[eta][phi][depth]=0;
00152                   subpedsum[eta][phi][depth]=0;
00153                   subpedsum2[eta][phi][depth]=0;
00154                   fC_rawpedsum[eta][phi][depth]=0;
00155                   fC_rawpedsum2[eta][phi][depth]=0;
00156                   fC_subpedsum[eta][phi][depth]=0;
00157                   fC_subpedsum2[eta][phi][depth]=0;
00158                 
00159                 } // loop over depths
00160             } // loop over phi
00161         } // loop over eta
00162 
00163     } // if (m_dbe)
00164 
00165   
00166   if (showTiming)
00167     {
00168       cpu_timer.stop();  cout <<"TIMER:: HcalPedestalMonitor SETUP -> "<<cpu_timer.cpuTime()<<endl;
00169     }
00170   
00171   return;
00172 
00173 } // void HcalPedestalMonitor::setup(...)
00174 
00175 
00176 // *************************************************** //
00177 
00178 
00179 void HcalPedestalMonitor::processEvent(const HBHEDigiCollection& hbhe,
00180                                        const HODigiCollection& ho,
00181                                        const HFDigiCollection& hf,
00182                                        //const ZDCDigiCollection& zdc, // ZDCs not yet added
00183                                        const HcalDbService& cond)
00184 {
00185   if (showTiming)
00186     {
00187       cpu_timer.reset(); cpu_timer.start();
00188     }
00189   
00190   ievt_++;
00191   meEVT_->Fill(ievt_);
00192 
00193   if(!shape_) shape_ = cond.getHcalShape(); // this one is generic
00194 
00195   if(!m_dbe) { 
00196     if(fVerbosity) cout<<"HcalPedestalMonitor::processEvent   DQMStore not instantiated!!!"<<endl;
00197     return; 
00198   }
00199   
00200   CaloSamples tool;  // digi values in ADC will be converted to fC values stored in tool
00201   float fC_myval=0;
00202   float ADC_myval=0;
00203 
00204   // HB/HE Loop
00205   try
00206     {    
00207       for (HBHEDigiCollection::const_iterator j=hbhe.begin(); 
00208            j!=hbhe.end(); ++j)
00209         {
00210       
00211           const HBHEDataFrame digi = (const HBHEDataFrame)(*j);
00212           if(!checkHB_ && (HcalSubdetector)(digi.id().subdet())==HcalBarrel) continue;
00213           if(!checkHE_ && (HcalSubdetector)(digi.id().subdet())==HcalEndcap) continue;
00214 
00215           calibs_= cond.getHcalCalibrations(digi.id());  // Old method was made private.
00216           int iEta   = digi.id().ieta();
00217           int iPhi   = digi.id().iphi();
00218           int iDepth = digi.id().depth();
00219 
00220           // Offset HBHE depths 1 and 2 by +4 so they don't overlap with HF
00221           if (digi.id().subdet()==HcalEndcap && iDepth!=3)
00222             iDepth+=4;
00223         
00224           channelCoder_ = cond.getHcalCoder(digi.id());
00225           HcalCoderDb coderDB(*channelCoder_, *shape_);
00226           coderDB.adc2fC(digi,tool);  // convert digi ADC counts to fC in tool
00227                   
00228           // Now loop over digi slices
00229           for (int k=0;k<digi.size();++k)
00230             {
00231               if (k<startingTimeSlice_ || k>endingTimeSlice_)
00232                 continue;
00233               
00234               unsigned int capid=digi.sample(k).capid();
00235               // Add ADC value to pedestal computation
00236               pedcounts[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]++;
00237               ADC_myval=digi.sample(k).adc();
00238               rawpedsum[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]+= ADC_myval;
00239               rawpedsum2[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]+=ADC_myval*ADC_myval;
00240               fC_rawpedsum[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]+= tool[k];
00241               fC_rawpedsum2[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]+=tool[k]*tool[k];
00242                       
00243               // Form subtracted pedestals
00244               // calibration object ALWAYS returns pedestals in fC, according to Radek (11 Feb 2009)
00245 
00246               fC_myval=tool[k]-calibs_.pedestal(capid);  // digi value in fC - pedestal (fC) for this capid
00247               //  HcalQIECoder->adc takes (shape, charge, capid) and returns adc value
00248               ADC_myval=digi.sample(k).adc()-(int)(channelCoder_->adc(*shape_,(float)calibs_.pedestal(capid), capid));
00249               
00250               //cout <<digi.sample(k).adc()<<"  "<<channelCoder_->adc(*shape_,(float)calibs_.pedestal(capid), capid)<<"  ADC = "<<ADC_myval<<endl;
00251               subpedsum[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]+=ADC_myval;
00252               subpedsum2[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]+=ADC_myval*ADC_myval;
00253               fC_subpedsum[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]+=fC_myval;
00254               fC_subpedsum2[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]+=fC_myval*fC_myval;
00255             } // for (int k=0;k<digi.size();++k)
00256         } // loop over digis
00257     } // try loop
00258   catch (...)
00259     {
00260       if (fVerbosity>0)
00261         cout <<"<HcalPedestalMonitor::processEvent>  No HBHE Digis."<<endl;
00262     }
00263   
00264   
00265   // HO Loop
00266 
00267   try
00268     {    
00269       for (HODigiCollection::const_iterator j=ho.begin(); 
00270            j!=ho.end(); ++j)
00271         {
00272           if (!checkHO_) continue;
00273           const HODataFrame digi = (const HODataFrame)(*j);
00274           //const HcalPedestalWidth* pedw = cond.getPedestalWidth(digi.id());
00275           calibs_= cond.getHcalCalibrations(digi.id());  // Old method was made private.
00276           int iEta   = digi.id().ieta();
00277           int iPhi   = digi.id().iphi();
00278           int iDepth = digi.id().depth();
00279           
00280           // Convert digi ADC value to fC (stored in tool)
00281           channelCoder_ = cond.getHcalCoder(digi.id());
00282           HcalCoderDb coderDB(*channelCoder_, *shape_);
00283           coderDB.adc2fC(digi,tool);
00284                       
00285           // Now loop over digi slices
00286           for (int k=0;k<digi.size();++k)
00287             {
00288               if (k<startingTimeSlice_ || k>endingTimeSlice_)
00289                 continue;
00290                   
00291               unsigned int capid=digi.sample(k).capid();
00292                   
00293               // Add ADC value to pedestal computation
00294               pedcounts[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]++;
00295               ADC_myval=digi.sample(k).adc();
00296               rawpedsum[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]+= ADC_myval;
00297               rawpedsum2[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]+=ADC_myval*ADC_myval;
00298               fC_rawpedsum[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]+= tool[k];
00299               fC_rawpedsum2[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]+=tool[k]*tool[k];
00300 
00301               // Form subtracted pedestals
00302               fC_myval=tool[k]-calibs_.pedestal(capid);
00303               //  HcalQIECoder->adc takes (shape, charge, capid) and returns adc value
00304               ADC_myval=digi.sample(k).adc()-(int)channelCoder_->adc(*shape_,(float)calibs_.pedestal(capid), capid);
00305 
00306               subpedsum[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]+=ADC_myval;
00307               subpedsum2[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]+=ADC_myval*ADC_myval;
00308               fC_subpedsum[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]+=fC_myval;
00309               fC_subpedsum2[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]+=fC_myval*fC_myval;
00310                   
00311             } // for (int k=0;k<digi.size();++k)
00312         } // loop over digis
00313     } // try loop
00314   catch (...)
00315     {
00316       if (fVerbosity > 0)
00317         cout <<"<HcalPedestalMonitor::processEvent>  No HO Digis."<<endl;
00318     }
00319 
00320 
00321   // HF Loop
00322   try
00323     {    
00324       for (HFDigiCollection::const_iterator j=hf.begin(); 
00325            j!=hf.end(); ++j)
00326         {
00327           if (!checkHO_) continue;
00328           const HFDataFrame digi = (const HFDataFrame)(*j);
00329           //const HcalPedestalWidth* pedw = cond.getPedestalWidth(digi.id());
00330           calibs_= cond.getHcalCalibrations(digi.id());  // Old method was made private.
00331           int iEta   = digi.id().ieta();
00332           int iPhi   = digi.id().iphi();
00333           int iDepth = digi.id().depth();
00334               
00335           channelCoder_ = cond.getHcalCoder(digi.id());
00336           HcalCoderDb coderDB(*channelCoder_, *shape_);
00337           coderDB.adc2fC(digi,tool);
00338                       
00339           // Now loop over digi slices
00340           for (int k=0;k<digi.size();++k)
00341             {
00342               if (k<startingTimeSlice_ || k>endingTimeSlice_)
00343                 continue;
00344 
00345               unsigned int capid=digi.sample(k).capid();
00346 
00347               // Add ADC value to pedestal computation
00348 
00349               // Depth values increased by 2 to avoid overlap with HE at |eta|=29
00350 
00351               pedcounts[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]++;
00352               ADC_myval=digi.sample(k).adc();
00353               rawpedsum[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]+= ADC_myval;
00354               rawpedsum2[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]+=ADC_myval*ADC_myval;
00355               fC_rawpedsum[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]+= tool[k];
00356               fC_rawpedsum2[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]+=tool[k]*tool[k];
00357 
00358               // calibs_.pedestal values are always in fC
00359               fC_myval=tool[k]-calibs_.pedestal(capid);
00360               //  HcalQIECoder->adc takes (shape, charge, capid) and returns adc value
00361               ADC_myval=digi.sample(k).adc()-(int)channelCoder_->adc(*shape_,(float)calibs_.pedestal(capid), capid);
00362 
00363               subpedsum[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]+=ADC_myval;
00364               subpedsum2[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]+=ADC_myval*ADC_myval;
00365               fC_subpedsum[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]+=fC_myval;
00366               fC_subpedsum2[iEta+(int)((etaBins_-2)/2)][iPhi-1][iDepth-1]+=fC_myval*fC_myval;
00367 
00368             } // for (int k=0;k<digi.size();++k)
00369         } // loop over digis
00370     } // try loop
00371   catch (...)
00372     {
00373       if (fVerbosity>0)
00374         cout <<"<HcalPedestalMonitor::processEvent>  No HF Digis."<<endl;
00375     }
00376 
00377   // Should we allow each subdetector to get filled separately?  (Fill hfHists every 1000 events, hoHists every 2000, etc.?)
00378 
00379   if (ievt_%pedmon_checkNevents_==0)
00380     {
00381       fillPedestalHistos();
00382     }
00383 
00384   if (showTiming)
00385     {
00386       cpu_timer.stop();  cout <<"TIMER:: HcalPedestalMonitor DIGI PROCESSEVENT -> "<<cpu_timer.cpuTime()<<endl;
00387     }
00388 
00389   
00390   return;
00391 } // void HcalPestalMonitor::processEvent(...)
00392 
00393 
00394 // *********************************************************** //
00395 
00396 
00397 void HcalPedestalMonitor::fillPedestalHistos(void)
00398 {
00399   if (showTiming)
00400     {
00401       cpu_timer.reset(); cpu_timer.start();
00402     }
00403 
00404   // Fills pedestal histograms
00405   if (fVerbosity>0) 
00406     cout <<"<HcalPedestalMonitor::fillPedestalHistos> Entered fillPedestalHistos routine"<<endl;
00407   
00408   // Set value to be filled in problem histograms to be checkNevents (or remainder of ievt_/pedmon_checkNevents_)
00409   
00410   double fillvalue=0;
00411 
00412   int mydepth=0;
00413   double ADC_myval, ADC_RMS, ADC_sub_myval, ADC_sub_RMS;
00414   
00415   double fC_myval, fC_RMS, fC_sub_myval, fC_sub_RMS;
00416   for (int eta=0;eta<(etaBins_-2);++eta)
00417     {
00418       for (int phi=0;phi<72;++phi)
00419         {
00420           for (int depth=0;depth<6;++depth) // this is one unit less "true" depth (for indexing purposes) 
00421             {
00422               mydepth=depth;
00423               
00424               // Skip events that don't contain required number of events
00425               if (pedcounts[eta][phi][depth] < minEntriesPerPed_) continue;
00426               
00427               // fillvalue = fraction of events used for pedestal determination
00428               // a small fillvalue causes the problem plots to get filled with a smaller value than a large fillvalue
00429               fillvalue = 1.*pedcounts[eta][phi][depth]/((endingTimeSlice_-startingTimeSlice_+1)*ievt_);
00430 
00431               // Compute mean and RMS for raw and subtracted pedestals in units of fC and ADC (phew!)
00432 
00433               ADC_myval= 1.*rawpedsum[eta][phi][depth]/pedcounts[eta][phi][depth]; 
00434               ADC_RMS = 1.0*rawpedsum2[eta][phi][depth]/pedcounts[eta][phi][depth]-1.*ADC_myval*ADC_myval;
00435               ADC_RMS=pow(fabs(ADC_RMS),0.5);
00436               
00437               ADC_sub_myval = 1.*subpedsum[eta][phi][depth]/pedcounts[eta][phi][depth];
00438               ADC_sub_RMS   = 1.0*subpedsum2[eta][phi][depth]/pedcounts[eta][phi][depth]-1.*ADC_sub_myval*ADC_sub_myval;
00439               ADC_sub_RMS = pow(fabs(ADC_sub_RMS),0.5);
00440                   
00441               fC_myval=1.*fC_rawpedsum[eta][phi][depth]/pedcounts[eta][phi][depth];
00442               fC_RMS=1.0*fC_rawpedsum2[eta][phi][depth]/pedcounts[eta][phi][depth]-1.*fC_myval*fC_myval;
00443               fC_RMS=pow(fabs(fC_RMS),0.5);
00444               
00445               /*
00446                 if ((eta-int((etaBins_-2)/2))==29 && depth==1)
00447                 cout <<phi<<"  "<<fC_subpedsum[eta][phi][depth]<<"  "<<fC_subpedsum2[eta][phi][depth]<<endl;
00448               */
00449               fC_sub_myval = 1.0*fC_subpedsum[eta][phi][depth]/pedcounts[eta][phi][depth];
00450               fC_sub_RMS   = 1.0*fC_subpedsum2[eta][phi][depth]/pedcounts[eta][phi][depth]-1.*fC_sub_myval*fC_sub_myval;
00451               fC_sub_RMS = pow(fabs(fC_sub_RMS),0.5);
00452 
00453               // When setting Bin Content, bins start at count of 1, not 0.
00454               // Also, first bins around eta,phi are empty.
00455               // Thus, eta,phi must be shifted by +2 (+1 for bin count, +1 to ignore empty row)
00456 
00457               // raw pedestals for HF
00458               rawADCPedestalMean[mydepth]->setBinContent(eta+2,phi+2,ADC_myval);
00459               rawADCPedestalRMS[mydepth]->setBinContent(eta+2,phi+2,ADC_RMS);
00460               rawADCPedestalMean_1D[mydepth]->Fill(ADC_myval);
00461               rawADCPedestalRMS_1D[mydepth]->Fill(ADC_RMS);
00462               rawFCPedestalMean[mydepth]->setBinContent(eta+2,phi+2,fC_myval);
00463               rawFCPedestalMean_1D[mydepth]->Fill(fC_myval);
00464               rawFCPedestalRMS[mydepth]->setBinContent(eta+2,phi+2,fC_RMS);
00465               rawFCPedestalRMS_1D[mydepth]->Fill(fC_RMS);
00466               
00467               //subtracted pedestals
00468               subADCPedestalMean[mydepth]->setBinContent(eta+2,phi+2,ADC_sub_myval);
00469               subADCPedestalRMS[mydepth]->setBinContent(eta+2,phi+2,ADC_sub_RMS);
00470               subADCPedestalMean_1D[mydepth]->Fill(ADC_sub_myval);
00471               subADCPedestalRMS_1D[mydepth]->Fill(ADC_sub_RMS);
00472               subFCPedestalMean[mydepth]->setBinContent(eta+2,phi+2,fC_sub_myval);
00473               subFCPedestalRMS[mydepth]->setBinContent(eta+2,phi+2,fC_sub_RMS);
00474               subFCPedestalMean_1D[mydepth]->Fill(fC_sub_myval);
00475               subFCPedestalRMS_1D[mydepth]->Fill(fC_sub_RMS);
00476               
00477               // Overall plots by depth
00478               MeanMapByDepth[mydepth]->setBinContent(eta+2,phi+2,ADC_myval);
00479               RMSMapByDepth[mydepth]->setBinContent(eta+2,phi+2,ADC_RMS);
00480               
00481               // Problem Cells
00482               // Problem Cells currently determined by checking ADC counts against
00483               // nominal expectations of (mean=3 ADC counts, RMS = 1 count)
00484               // We may instead want to compare the values to the database values in order to determine problems?
00485               if  (fillvalue>pedmon_minErrorFlag_ 
00486                    && (fabs(ADC_myval-nominalPedMeanInADC_)>maxPedMeanDiffADC_ 
00487                        || fabs(ADC_RMS-nominalPedWidthInADC_)>maxPedWidthDiffADC_))
00488                 {
00489                   ProblemPedestals->setBinContent(eta+2,phi+2,fillvalue);
00490                   ProblemPedestalsByDepth[mydepth]->setBinContent(eta+2,phi+2,fillvalue);
00491                 }
00492               
00493                   // ZDC still to be added
00494                   /*
00495                     ADD ZDC HERE AT SOME POINT!
00496                   */
00497               
00498             } // for (int depth)
00499         } // for (int phi)
00500     } // for (int eta)
00501   
00502   // Fill unphysical cells
00503   FillUnphysicalHEHFBins(rawADCPedestalMean);
00504   FillUnphysicalHEHFBins(rawADCPedestalRMS);
00505   FillUnphysicalHEHFBins(rawFCPedestalMean); 
00506   FillUnphysicalHEHFBins(rawFCPedestalRMS);  
00507   
00508   //subtracted pedestals
00509   FillUnphysicalHEHFBins(subADCPedestalMean);
00510   FillUnphysicalHEHFBins(subADCPedestalRMS); 
00511   
00512   FillUnphysicalHEHFBins(subFCPedestalMean); 
00513   FillUnphysicalHEHFBins(subFCPedestalRMS);  
00514   
00515   // Overall plots by depth
00516   FillUnphysicalHEHFBins(MeanMapByDepth);    
00517   FillUnphysicalHEHFBins(RMSMapByDepth);     
00518   
00519   FillUnphysicalHEHFBins(ProblemPedestalsByDepth); 
00520   FillUnphysicalHEHFBins(ProblemPedestals);
00521 
00522   if (showTiming)
00523     {
00524       cpu_timer.stop();  cout <<"TIMER:: HcalPedestalHistos DIGI FILLPEDESTALHISTOS -> "<<cpu_timer.cpuTime()<<endl;
00525     }
00526 
00527   return;
00528 
00529 } //void HcalPedestalMonitor::fillPedestalHistos(void)
00530 
00531 
00532 // *********************************************************** //
00533 
00534 
00535 void HcalPedestalMonitor::done()
00536 {
00537   // I'd like to put in another call to fillPedestalHistos() here, but this gets called after root file gets written?
00538   // update on 22 October 2008:  DQMFileSaver gets called with endRun, not endJob.  If we want the the plots to be updated, we should call them in endRun
00539   return;
00540 }
00541 
00542 
00543 // ******************************************************** //
00544 
00545 void HcalPedestalMonitor::fillDBValues(const HcalDbService& cond)
00546 {
00547   /* Fills reference pedestal mean, width plots with pedestal values from conditions database
00548    */
00549   
00550   if (showTiming)
00551     {
00552       cpu_timer.reset(); cpu_timer.start();
00553     }
00554 
00555 
00556   if(!shape_) shape_ = cond.getHcalShape(); // this one is generic
00557   
00558   int fill_offset=0;
00559   double ADC_ped=0;
00560   double ADC_width=0;
00561   double fC_ped=0;
00562   double fC_width=0;
00563   for (int subdet=1; subdet<=4;++subdet)
00564     {
00565       for (int depth=1;depth<=4;++depth)
00566         {
00567           for (int ieta=(int)etaMin_;ieta<=(int)etaMax_;++ieta)
00568             {
00569               for (int iphi=(int)phiMin_;iphi<=(int)phiMax_;++iphi)
00570                 {
00571                   //if (!hcal.validDetId((HcalSubdetector)(subdet), ieta, iphi, depth)) continue; // implement once this is available in future version of HcalDetId.h
00572                   if (!validDetId((HcalSubdetector)(subdet), ieta, iphi, depth)) continue;
00573                   HcalDetId detid((HcalSubdetector)(subdet), ieta, iphi, depth);
00574 
00575                   ADC_ped=0;
00576                   ADC_width=0;
00577                   fC_ped=0;
00578                   fC_width=0;
00579                   calibs_= cond.getHcalCalibrations(detid);  
00580                   const HcalPedestalWidth* pedw = cond.getPedestalWidth(detid);
00581                   channelCoder_ = cond.getHcalCoder(detid);
00582                   
00583                   // Loop over capIDs
00584                   for (unsigned int capid=0;capid<4;++capid)
00585                     {
00586                       // Still need to determine how to convert widths to ADC or fC
00587                       // calibs_.pedestal value is always in fC, according to Radek
00588                       fC_ped+=calibs_.pedestal(capid);
00589                       // convert to ADC from fC
00590                       ADC_ped+=channelCoder_->adc(*shape_,
00591                                                   (float)calibs_.pedestal(capid),
00592                                                   capid);
00593 
00594                       if (doFCpeds_)
00595                         {
00596                           fC_width+=pedw->getSigma(capid,capid);
00597                           ADC_width+=pedw->getSigma(capid,capid)*pow(1.*channelCoder_->adc(*shape_,(float)calibs_.pedestal(capid),capid)/calibs_.pedestal(capid),2);
00598                         }
00599                       else
00600                         {
00601                           ADC_width+=pedw->getSigma(capid,capid);
00602                           fC_width+=pedw->getSigma(capid,capid)*pow(1.*calibs_.pedestal(capid)/channelCoder_->adc(*shape_,(float)calibs_.pedestal(capid),capid),2);
00603                         }
00604                     }//capid loop
00605 
00606                   // Pedestal values are average over four cap IDs
00607                   // widths are sqrt(SUM [sigma_ii^2])/4.
00608                   fC_ped/=4.;
00609                   ADC_ped/=4.;
00610                   fC_width=pow(fC_width,0.5)/4.;
00611                   ADC_width=pow(ADC_width,0.5)/4.;
00612 
00613                   if (fVerbosity>1)
00614                     {
00615                       cout <<"<HcalPedestalMonitor::fillDBValues> HcalDet ID = "<<(HcalSubdetector)subdet<<": ("<<ieta<<", "<<iphi<<", "<<depth<<")"<<endl;
00616                       cout <<"\tADC pedestal = "<<ADC_ped<<" +/- "<<ADC_width<<endl;
00617                       cout <<"\tfC pedestal = "<<fC_ped<<" +/- "<<fC_width<<endl;
00618                     }
00619                   ADC_PedestalFromDB->Fill(ieta,iphi,ADC_ped);
00620                   ADC_WidthFromDB->Fill(ieta,iphi,ADC_width);
00621                   fC_PedestalFromDB->Fill(ieta,iphi,fC_ped);
00622                   fC_WidthFromDB->Fill(ieta,iphi,fC_width);
00623 
00624                   if (((HcalSubdetector)(subdet)==HcalEndcap) && depth<3) fill_offset=4;
00625                   else fill_offset=0;
00626                   ADC_PedestalFromDBByDepth[depth-1+fill_offset]->Fill(ieta,iphi,ADC_ped);
00627                   ADC_WidthFromDBByDepth[depth-1+fill_offset]->Fill(ieta, iphi, ADC_width);
00628                   fC_PedestalFromDBByDepth[depth-1+fill_offset]->Fill(ieta,iphi,fC_ped);
00629                   fC_WidthFromDBByDepth[depth-1+fill_offset]->Fill(ieta, iphi, fC_width);
00630                 } // iphi loop
00631             } // ieta loop
00632         } //depth loop
00633 
00634     } // subdet loop
00635    
00636   if (showTiming)
00637     {
00638       cpu_timer.stop();  
00639       cout <<"TIMER:: HcalPedestalMonitor FILLDBVALUES -> "<<cpu_timer.cpuTime()<<endl;
00640     }
00641 
00642   return;
00643 } // void HcalPedestalMonitor::fillDBValues(void)
00644 
00645 
00646 // ******************************************************** //

Generated on Tue Jun 9 17:33:02 2009 for CMSSW by  doxygen 1.5.4