CMS 3D CMS Logo

HcalMonitorModule.cc

Go to the documentation of this file.
00001 #include <DQM/HcalMonitorModule/src/HcalMonitorModule.h>
00002 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
00003 
00004 /*
00005  * \file HcalMonitorModule.cc
00006  * 
00007  * $Date: 2009/04/05 13:31:15 $
00008  * $Revision: 1.111.2.1 $
00009  * \author W Fisher
00010  * \author J Temple
00011  *
00012 */
00013 
00014 using namespace std;
00015 using namespace edm;
00016 
00017 //--------------------------------------------------------
00018 HcalMonitorModule::HcalMonitorModule(const edm::ParameterSet& ps){
00019 
00020   irun_=0; ilumisec_=0; ievent_=0; itime_=0;
00021   actonLS_=false;
00022   meStatus_=0;  meRunType_=0;
00023   meEvtMask_=0; meFEDS_=0;
00024   meLatency_=0; meQuality_=0;
00025   fedsListed_ = false;
00026   digiMon_ = 0;   dfMon_ = 0;
00027   diTask_ = 0;
00028   rhMon_ = 0;     pedMon_ = 0; 
00029   ledMon_ = 0;    mtccMon_ = 0;
00030   hotMon_ = 0;    tempAnalysis_ = 0;
00031   deadMon_ = 0;   tpMon_ = 0;
00032   ctMon_ = 0;     beamMon_ = 0;
00033   laserMon_ = 0;
00034   expertMon_ = 0;  eeusMon_ = 0;
00035 
00036   // initialize hcal quality object
00037   
00038 
00039   // All subdetectors assumed out of the run by default
00040   HBpresent_=0;
00041   HEpresent_=0;
00042   HOpresent_=0;
00043   HFpresent_=0;
00044 
00045   inputLabelDigi_        = ps.getParameter<edm::InputTag>("digiLabel");
00046   inputLabelRecHitHBHE_  = ps.getParameter<edm::InputTag>("hbheRecHitLabel");
00047   inputLabelRecHitHF_    = ps.getParameter<edm::InputTag>("hfRecHitLabel");
00048   inputLabelRecHitHO_    = ps.getParameter<edm::InputTag>("hoRecHitLabel");
00049   inputLabelRecHitZDC_   = ps.getParameter<edm::InputTag>("zdcRecHitLabel");
00050   inputLabelCaloTower_   = ps.getParameter<edm::InputTag>("caloTowerLabel");
00051   inputLabelLaser_       = ps.getParameter<edm::InputTag>("hcalLaserLabel");
00052 
00053   checkHB_=ps.getUntrackedParameter<bool>("checkHB", 1); 
00054   checkHE_=ps.getUntrackedParameter<bool>("checkHE", 1);  
00055   checkHO_=ps.getUntrackedParameter<bool>("checkHO", 1);  
00056   checkHF_=ps.getUntrackedParameter<bool>("checkHF", 1);   
00057 
00058   evtSel_ = new HcalMonitorSelector(ps);
00059   
00060   dbe_ = Service<DQMStore>().operator->();
00061   
00062   debug_ = ps.getUntrackedParameter<int>("debug", 0);
00063   
00064   showTiming_ = ps.getUntrackedParameter<bool>("showTiming", false);
00065   dump2database_   = ps.getUntrackedParameter<bool>("dump2database",false); // dumps output to database file
00066 
00067   FEDRawDataCollection_ = ps.getUntrackedParameter<edm::InputTag>("FEDRawDataCollection",edm::InputTag("source",""));
00068 
00069   // Valgrind complained when the test was simply:  if ( ps.getUntrackedParameter<bool>("DataFormatMonitor", false))
00070   // try assigning value to bool first?
00071   bool taskOn = ps.getUntrackedParameter<bool>("DataFormatMonitor", false);
00072   if (taskOn) {
00073     if(debug_>0) cout << "HcalMonitorModule: DataFormat monitor flag is on...." << endl;
00074     dfMon_ = new HcalDataFormatMonitor();
00075     dfMon_->setup(ps, dbe_);
00076   }
00077 
00078   taskOn = ps.getUntrackedParameter<bool>("DataIntegrityTask", false); 
00079   if (taskOn ) 
00080     {
00081       if (debug_>0) cout <<"HcalMonitorModule: DataIntegrity monitor flag is on...."<<endl;
00082       diTask_ = new HcalDataIntegrityTask();
00083       diTask_->setup(ps, dbe_);
00084     }
00085 
00086   if ( ps.getUntrackedParameter<bool>("DigiMonitor", false) ) {
00087     if(debug_>0) cout << "HcalMonitorModule: Digi monitor flag is on...." << endl;
00088     digiMon_ = new HcalDigiMonitor();
00089     digiMon_->setup(ps, dbe_);
00090   }
00091   
00092   if ( ps.getUntrackedParameter<bool>("RecHitMonitor", false) ) {
00093     if(debug_>0) cout << "HcalMonitorModule: RecHit monitor flag is on...." << endl;
00094     rhMon_ = new HcalRecHitMonitor();
00095     rhMon_->setup(ps, dbe_);
00096   }
00097   
00098   if ( ps.getUntrackedParameter<bool>("PedestalMonitor", false) ) {
00099     if(debug_>0) cout << "HcalMonitorModule: Pedestal monitor flag is on...." << endl;
00100     pedMon_ = new HcalPedestalMonitor();
00101     pedMon_->setup(ps, dbe_);
00102   }
00103   
00104   if ( ps.getUntrackedParameter<bool>("LEDMonitor", false) ) {
00105     if(debug_>0) cout << "HcalMonitorModule: LED monitor flag is on...." << endl;
00106     ledMon_ = new HcalLEDMonitor();
00107     ledMon_->setup(ps, dbe_);
00108   }
00109   
00110   if ( ps.getUntrackedParameter<bool>("LaserMonitor", false) ) {
00111     if(debug_>0) cout << "HcalMonitorModule: Laser monitor flag is on...." << endl;
00112     laserMon_ = new HcalLaserMonitor();
00113     laserMon_->setup(ps, dbe_);
00114   }
00115 
00116   if ( ps.getUntrackedParameter<bool>("MTCCMonitor", false) ) {
00117     if(debug_>0) cout << "HcalMonitorModule: MTCC monitor flag is on...." << endl;
00118     mtccMon_ = new HcalMTCCMonitor();
00119     mtccMon_->setup(ps, dbe_);
00120   }
00121   
00122   if ( ps.getUntrackedParameter<bool>("HotCellMonitor", false) ) {
00123     if(debug_>0) cout << "HcalMonitorModule: Hot Cell monitor flag is on...." << endl;
00124     hotMon_ = new HcalHotCellMonitor();
00125     hotMon_->setup(ps, dbe_);
00126   }
00127   
00128   if ( ps.getUntrackedParameter<bool>("DeadCellMonitor", false) ) {
00129     if(debug_>0) cout << "HcalMonitorModule: Dead Cell monitor flag is on...." << endl;
00130     deadMon_ = new HcalDeadCellMonitor();
00131     deadMon_->setup(ps, dbe_);
00132   }
00133 
00134   if ( ps.getUntrackedParameter<bool>("TrigPrimMonitor", false) ) {      
00135     if(debug_>0) cout << "HcalMonitorModule: TrigPrim monitor flag is on...." << endl;   
00136     tpMon_ = new HcalTrigPrimMonitor();          
00137     tpMon_->setup(ps, dbe_);     
00138   }  
00139 
00140   if (ps.getUntrackedParameter<bool>("CaloTowerMonitor",false)){
00141     if(debug_>0) cout << "HcalMonitorModule: CaloTower monitor flag is on...." << endl;          
00142     ctMon_ = new HcalCaloTowerMonitor();         
00143     ctMon_->setup(ps, dbe_);     
00144   }  
00145 
00146   if (ps.getUntrackedParameter<bool>("BeamMonitor",false)){
00147     if(debug_>0) cout << "HcalMonitorModule: Beam monitor flag is on...."<<endl;
00148     beamMon_ = new HcalBeamMonitor();
00149     beamMon_->setup(ps, dbe_);
00150   }
00151 
00152   if (ps.getUntrackedParameter<bool>("ExpertMonitor",false)){
00153     if(debug_>0) cout << "HcalMonitorModule: Expert monitor flag is on...."<<endl;
00154     expertMon_ = new HcalExpertMonitor();
00155     expertMon_->setup(ps, dbe_);
00156   }
00157 
00158   
00159 
00160   if ( ps.getUntrackedParameter<bool>("HcalAnalysis", false) ) {
00161     if(debug_>0) cout << "HcalMonitorModule: Hcal Analysis flag is on...." << endl;
00162     tempAnalysis_ = new HcalTemplateAnalysis();
00163     tempAnalysis_->setup(ps);
00164   }
00165 
00166   if (ps.getUntrackedParameter<bool>("EEUSMonitor",false))
00167     {
00168       if (debug_>0) cout <<"HcalMonitorModule:  Empty Event/Unsuppressed Moniotr is on..."<<endl;
00169       eeusMon_ = new HcalEEUSMonitor();
00170       eeusMon_->setup(ps, dbe_);
00171     }
00172 
00173   // set parameters   
00174   prescaleEvt_ = ps.getUntrackedParameter<int>("diagnosticPrescaleEvt", -1);
00175   if(debug_>1) cout << "===>HcalMonitor event prescale = " << prescaleEvt_ << " event(s)"<< endl;
00176 
00177   prescaleLS_ = ps.getUntrackedParameter<int>("diagnosticPrescaleLS", -1);
00178   if(debug_>1) cout << "===>HcalMonitor lumi section prescale = " << prescaleLS_ << " lumi section(s)"<< endl;
00179   if (prescaleLS_>0) actonLS_=true;
00180 
00181   prescaleUpdate_ = ps.getUntrackedParameter<int>("diagnosticPrescaleUpdate", -1);
00182   if(debug_>1) cout << "===>HcalMonitor update prescale = " << prescaleUpdate_ << " update(s)"<< endl;
00183 
00184   prescaleTime_ = ps.getUntrackedParameter<int>("diagnosticPrescaleTime", -1);
00185   if(debug_>1) cout << "===>HcalMonitor time prescale = " << prescaleTime_ << " minute(s)"<< endl;
00186   
00187   // Base folder for the contents of this job
00188   string subsystemname = ps.getUntrackedParameter<string>("subSystemFolder", "Hcal") ;
00189   if(debug_>0) cout << "===>HcalMonitor name = " << subsystemname << endl;
00190   rootFolder_ = subsystemname + "/";
00191   
00192   gettimeofday(&psTime_.updateTV,NULL);
00194   psTime_.updateTime = (psTime_.updateTV.tv_sec*1000.0+psTime_.updateTV.tv_usec/1000.0);
00195   psTime_.updateTime /= 1000.0;
00196   psTime_.elapsedTime=0;
00197   psTime_.vetoTime=psTime_.updateTime;
00198 }
00199 
00200 //--------------------------------------------------------
00201 HcalMonitorModule::~HcalMonitorModule()
00202 {
00203   
00204   if (dbe_!=0)
00205     {    
00206       if(digiMon_!=0)   {  digiMon_->clearME();}
00207      if(dfMon_!=0)     {  dfMon_->clearME();}
00208      if(diTask_!=0)    {  diTask_->clearME();}
00209      if(pedMon_!=0)    {  pedMon_->clearME();}
00210      if(ledMon_!=0)    {  ledMon_->clearME();}
00211      if(laserMon_!=0)  {  laserMon_->clearME();}
00212      if(hotMon_!=0)    {  hotMon_->clearME();}
00213      if(deadMon_!=0)   {  deadMon_->clearME();}
00214      if(mtccMon_!=0)   {  mtccMon_->clearME();}
00215      if(rhMon_!=0)     {  rhMon_->clearME();}
00216      
00217      dbe_->setCurrentFolder(rootFolder_);
00218      dbe_->removeContents();
00219     }
00220   
00221   // I think setting pointers to NULL (0) after delete is unnecessary here,
00222   // since we're in the destructor (and thus won't be using the pointers again.)
00223   if(digiMon_!=0) 
00224     { 
00225       delete digiMon_;  digiMon_=0; 
00226     }
00227   if(dfMon_!=0) 
00228     { delete dfMon_;     dfMon_=0; 
00229     }
00230   if(diTask_!=0) 
00231     { delete diTask_;   diTask_=0; 
00232     }
00233   if(pedMon_!=0) 
00234     {
00235       delete pedMon_;   pedMon_=0; 
00236     }
00237   if(ledMon_!=0) 
00238     { delete ledMon_;   ledMon_=0; 
00239     }
00240   if(laserMon_!=0) 
00241     { delete laserMon_;   laserMon_=0; 
00242     }
00243   if(hotMon_!=0) 
00244     { delete hotMon_;   hotMon_=0; 
00245     }
00246   if(deadMon_!=0) 
00247     { delete deadMon_; deadMon_=0; 
00248     }
00249   if (beamMon_!=0)
00250     { delete beamMon_;  beamMon_=0;
00251     }
00252   if(mtccMon_!=0) 
00253     { delete mtccMon_; mtccMon_=0; 
00254     }
00255   if(rhMon_!=0) 
00256     { delete rhMon_;     rhMon_=0; 
00257     }
00258   if(tempAnalysis_!=0) 
00259     { delete tempAnalysis_; 
00260     tempAnalysis_=0; 
00261     }
00262   if (evtSel_!=0) {delete evtSel_; evtSel_ = 0;
00263   }
00264 } //void HcalMonitorModule::~HcalMonitorModule()
00265 
00266 //--------------------------------------------------------
00267 void HcalMonitorModule::beginJob(const edm::EventSetup& c){
00268   ievt_ = 0;
00269   
00270   ievt_pre_=0;
00271 
00272   if ( dbe_ != NULL ){
00273     dbe_->setCurrentFolder(rootFolder_+"DQM Job Status" );
00274     meStatus_  = dbe_->bookInt("STATUS");
00275     meRunType_ = dbe_->bookInt("RUN TYPE");
00276     meEvtMask_ = dbe_->bookInt("EVT MASK");
00277     meFEDS_    = dbe_->book1D("FEDs Unpacked","FEDs Unpacked",100,700,799);
00278     // process latency was (200,0,1), but that gave overflows
00279     meLatency_ = dbe_->book1D("Process Latency","Process Latency",2000,0,10);
00280     meQuality_ = dbe_->book1D("Quality Status","Quality Status",100,0,1);
00281     // Store whether or not subdetectors are present
00282     meHB_ = dbe_->bookInt("HBpresent");
00283     meHE_ = dbe_->bookInt("HEpresent");
00284     meHO_ = dbe_->bookInt("HOpresent");
00285     meHF_ = dbe_->bookInt("HFpresent");
00286     
00287     meStatus_->Fill(0);
00288     meRunType_->Fill(-1);
00289     meEvtMask_->Fill(-1);
00290 
00291     // Should fill with 0 to start
00292     meHB_->Fill(HBpresent_);
00293     meHE_->Fill(HEpresent_);
00294     meHO_->Fill(HOpresent_);
00295     meHF_->Fill(HFpresent_);
00296   }
00297 
00298   edm::ESHandle<HcalDbService> pSetup;
00299   c.get<HcalDbRecord>().get( pSetup );
00300 
00301   readoutMap_=pSetup->getHcalMapping();
00302   DetId detid_;
00303   HcalDetId hcaldetid_; 
00304 
00305   // Build a map of readout hardware unit to calorimeter channel
00306   std::vector <HcalElectronicsId> AllElIds = readoutMap_->allElectronicsIdPrecision();
00307   int dccid;
00308   pair <int,int> dcc_spgt;
00309   // by looping over all precision (non-trigger) items.
00310   for (std::vector <HcalElectronicsId>::iterator eid = AllElIds.begin();
00311        eid != AllElIds.end();
00312        eid++) {
00313 
00314     //Get the HcalDetId from the HcalElectronicsId
00315     detid_ = readoutMap_->lookup(*eid);
00316     
00317 
00318     // NULL if illegal; ignore
00319     if (!detid_.null()) {
00320       if (detid_.det()!=4) continue;
00321       if (detid_.subdetId()!=HcalBarrel &&
00322           detid_.subdetId()!=HcalEndcap &&
00323           detid_.subdetId()!=HcalOuter  &&
00324           detid_.subdetId()!=HcalForward) continue;
00325       hcaldetid_ = HcalDetId(detid_);
00326       
00327       dccid = eid->dccid();
00328       dcc_spgt = pair <int,int> (dccid, eid->spigot());
00329       
00330       thisDCC = DCCtoCell.find(dccid);
00331       thisHTR = HTRtoCell.find(dcc_spgt);
00332       
00333       // If this DCC has no entries, make this its first one.
00334       if (thisDCC == DCCtoCell.end()) {
00335         std::vector <HcalDetId> tempv;
00336         tempv.push_back(hcaldetid_);
00337         pair <int, std::vector<HcalDetId> > thispair;
00338         thispair = pair <int, std::vector<HcalDetId> > (dccid,tempv);
00339         DCCtoCell.insert(thispair); 
00340       }
00341       else {
00342         thisDCC->second.push_back(hcaldetid_);
00343       }
00344       
00345       // If this HTR has no entries, make this its first one.
00346       if (thisHTR == HTRtoCell.end()) {
00347         std::vector <HcalDetId> tempv;
00348         tempv.push_back(hcaldetid_);
00349         pair < pair <int,int>, std::vector<HcalDetId> > thispair;
00350         thispair = pair <pair <int,int>, std::vector<HcalDetId> > (dcc_spgt,tempv);
00351         HTRtoCell.insert(thispair); 
00352       }
00353       else {
00354         thisHTR->second.push_back(hcaldetid_);  
00355       }
00356     } // if (!detid_.null()) 
00357   } 
00358   if (dfMon_) {
00359     dfMon_->smuggleMaps(DCCtoCell, HTRtoCell);
00360   }
00361 
00362   //get conditions
00363   c.get<HcalDbRecord>().get(conditions_);
00364 
00365   // fill reference pedestals with database values
00366   // Need to repeat this so many times?  Just do it once? And then we can be smarter about the whole fC/ADC thing?
00367   if (pedMon_!=NULL)
00368     pedMon_->fillDBValues(*conditions_);
00369   if (deadMon_!=NULL)
00370     deadMon_->createMaps(*conditions_);
00371   if (hotMon_!=NULL)
00372     hotMon_->createMaps(*conditions_);
00373 
00374 
00375   edm::ESHandle<HcalChannelQuality> p;
00376   c.get<HcalChannelQualityRcd>().get(p);
00377   chanquality_= new HcalChannelQuality(*p.product());
00378   return;
00379 } // HcalMonitorModule::beginJob(...)
00380 
00381 //--------------------------------------------------------
00382 void HcalMonitorModule::beginRun(const edm::Run& run, const edm::EventSetup& c) {
00383   fedsListed_ = false;
00384 
00385   // I think we want to reset these at 0 at the start of each run
00386   HBpresent_ = 0;
00387   HEpresent_ = 0;
00388   HOpresent_ = 0;
00389   HFpresent_ = 0;
00390 
00391   // Should fill with 0 to start
00392   meHB_->Fill(HBpresent_);
00393   meHE_->Fill(HEpresent_);
00394   meHO_->Fill(HOpresent_);
00395   meHF_->Fill(HFpresent_);
00396   reset();
00397 }
00398 
00399 //--------------------------------------------------------
00400 void HcalMonitorModule::beginLuminosityBlock(const edm::LuminosityBlock& lumiSeg, 
00401      const edm::EventSetup& context) {
00402   
00403   if(actonLS_ && !prescale()){
00404     // do scheduled tasks...
00405   }
00406 }
00407 
00408 
00409 //--------------------------------------------------------
00410 void HcalMonitorModule::endLuminosityBlock(const edm::LuminosityBlock& lumiSeg, 
00411                                            const edm::EventSetup& context) {
00412   if(actonLS_ && !prescale()){
00413     // do scheduled tasks...
00414   }
00415 }
00416 
00417 //--------------------------------------------------------
00418 void HcalMonitorModule::endRun(const edm::Run& r, const edm::EventSetup& context)
00419 {
00420   if (debug_>0)  
00421     cout <<"HcalMonitorModule::endRun(...) ievt = "<<ievt_<<endl;
00422 
00423   // Do final pedestal histogram filling
00424   if (pedMon_!=NULL)
00425     pedMon_->fillPedestalHistos();
00426 
00427   if (deadMon_!=NULL)
00428     deadMon_->fillDeadHistosAtEndRun();
00429 
00430   return;
00431 }
00432 
00433 
00434 //--------------------------------------------------------
00435 void HcalMonitorModule::endJob(void) {
00436   
00437   if ( meStatus_ ) meStatus_->Fill(2);
00438 
00439   if(rhMon_!=NULL) rhMon_->done();
00440   if(digiMon_!=NULL) digiMon_->done();
00441   if(dfMon_!=NULL) dfMon_->done();
00442   if(diTask_!=NULL) diTask_->done();
00443   if(pedMon_!=NULL) pedMon_->done();
00444   if(ledMon_!=NULL) ledMon_->done();
00445   if(laserMon_!=NULL) laserMon_->done();
00446   if(hotMon_!=NULL) hotMon_->done(myquality_);
00447   if(deadMon_!=NULL) deadMon_->done(myquality_);
00448   if(mtccMon_!=NULL) mtccMon_->done();
00449   if (tpMon_!=NULL) tpMon_->done();
00450   if (ctMon_!=NULL) ctMon_->done();
00451   if (beamMon_!=NULL) beamMon_->done();
00452   if (expertMon_!=NULL) expertMon_->done();
00453   if (eeusMon_!=NULL) eeusMon_->done();
00454   if(tempAnalysis_!=NULL) tempAnalysis_->done();
00455 
00456   if (dump2database_)
00457     {
00458       if (debug_>0) cout <<"<HcalMonitorModule::endJob>  Writing file for database"<<endl;
00459       std::vector<DetId> mydetids = chanquality_->getAllChannels();
00460       HcalChannelQuality* newChanQual = new HcalChannelQuality();
00461       for (unsigned int i=0;i<mydetids.size();++i)
00462         {
00463           if (mydetids[i].det()!=4) continue; // not hcal
00464           //HcalDetId id(mydetids[i]);
00465           HcalDetId id=mydetids[i];
00466           // get original channel status item
00467           const HcalChannelStatus* origstatus=chanquality_->getValues(mydetids[i]);
00468           // make copy of status
00469           HcalChannelStatus* mystatus=new HcalChannelStatus(origstatus->rawId(),origstatus->getValue());
00470           if (myquality_.find(id)!=myquality_.end())
00471             {
00472               // Set bit 1 for cells which aren't present        
00473               if ((id.subdet()==HcalBarrel &&!HBpresent_) ||     
00474                   (id.subdet()==HcalEndcap &&!HEpresent_) ||     
00475                   (id.subdet()==HcalOuter  &&!HOpresent_) ||     
00476                   (id.subdet()==HcalForward&&!HFpresent_))       
00477                 {        
00478                   mystatus->setBit(1);   
00479                 }        
00480               // Only perform these checks if bit 0 not set?
00481               // check dead cells
00482               if ((myquality_[id]>>5)&0x1)
00483                   mystatus->setBit(5);
00484               else
00485                 mystatus->unsetBit(5);
00486               // check hot cells
00487               if ((myquality_[id]>>6)&0x1)
00488                 mystatus->setBit(6);
00489               else
00490                 mystatus->unsetBit(6);
00491             } // if (myquality_.find_...)
00492           newChanQual->addValues(*mystatus);
00493           // Clean up pointers to avoid memory leaks
00494           delete origstatus;
00495           delete mystatus;
00496         } // for (unsigned int i=0;...)
00497       // Now dump out to text file
00498       std::ostringstream file;
00499       file <<"HcalDQMstatus_"<<irun_<<".txt";
00500       std::ofstream outStream(file.str().c_str());
00501       HcalDbASCIIIO::dumpObject (outStream, (*newChanQual));
00502 
00503     } // if (dump2databse_)
00504   return;
00505 }
00506 
00507 //--------------------------------------------------------
00508 void HcalMonitorModule::reset(){
00509 
00510   if(rhMon_!=NULL)   rhMon_->reset();
00511   if(digiMon_!=NULL) digiMon_->reset();
00512   if(dfMon_!=NULL)   dfMon_->reset();
00513   if(diTask_!=NULL)  diTask_->reset();
00514   if(pedMon_!=NULL)  pedMon_->reset();
00515   if(ledMon_!=NULL)  ledMon_->reset();
00516   if(laserMon_!=NULL)  laserMon_->reset();
00517   if(hotMon_!=NULL)  hotMon_->reset();
00518   if(deadMon_!=NULL)  deadMon_->reset();
00519   if(mtccMon_!=NULL)   mtccMon_->reset();
00520   if(tempAnalysis_!=NULL) tempAnalysis_->reset();
00521   if(tpMon_!=NULL) tpMon_->reset();
00522   if(ctMon_!=NULL) ctMon_->reset();
00523   if(beamMon_!=NULL) beamMon_->reset();
00524   if(expertMon_!=NULL) expertMon_->reset();
00525   if(eeusMon_!=NULL) eeusMon_->reset();
00526 
00527 }
00528 
00529 //--------------------------------------------------------
00530 void HcalMonitorModule::analyze(const edm::Event& e, const edm::EventSetup& eventSetup){
00531 
00532   // environment datamembers
00533   irun_     = e.id().run();
00534   ilumisec_ = e.luminosityBlock();
00535   ievent_   = e.id().event();
00536   itime_    = e.time().value();
00537 
00538   if (debug_>1) cout << "HcalMonitorModule: evts: "<< nevt_ << ", run: " << irun_ << ", LS: " << ilumisec_ << ", evt: " << ievent_ << ", time: " << itime_ << endl <<"\t counter = "<<ievt_pre_<<"\t total count = "<<ievt_<<endl; 
00539 
00540   // skip this event if we're prescaling...
00541   ievt_pre_++; // need to increment counter before calling prescale
00542 
00543   if(prescale()) return;
00544   meLatency_->Fill(psTime_.elapsedTime);
00545 
00546   // Do default setup...
00547   ievt_++;
00548 
00549   int evtMask=DO_HCAL_DIGIMON|DO_HCAL_DFMON|DO_HCAL_RECHITMON|DO_HCAL_PED_CALIBMON|DO_HCAL_LED_CALIBMON|DO_HCAL_LASER_CALIBMON; // add in DO_HCAL_TPMON, DO_HCAL_CTMON?  (in HcalMonitorSelector.h)
00550 
00551   //  int trigMask=0;
00552   if(mtccMon_==NULL){
00553     evtSel_->processEvent(e);
00554     evtMask = evtSel_->getEventMask();
00555     //    trigMask =  evtSel_->getTriggerMask();
00556   }
00557   if ( dbe_ ){ 
00558     meStatus_->Fill(1);
00559     meEvtMask_->Fill(evtMask);
00560   }
00561   
00563   bool rawOK_    = true;
00564   bool digiOK_   = true;
00565   bool rechitOK_ = true;
00566   bool zdchitOK_ = true;
00567   bool trigOK_   = false;
00568   bool tpdOK_    = true;
00569   bool calotowerOK_ = true;
00570   bool laserOK_  = true;
00571 
00572   // try to get raw data and unpacker report
00573   edm::Handle<FEDRawDataCollection> rawraw;  
00574 
00575   // Trying new getByLabel
00576   if (!(e.getByLabel(FEDRawDataCollection_,rawraw)))
00577     {
00578       rawOK_=false;
00579       LogWarning("HcalMonitorModule")<<" source not available";
00580     }
00581   if (rawOK_&&!rawraw.isValid()) {
00582     rawOK_=false;
00583   }
00584 
00585   edm::Handle<HcalUnpackerReport> report;  
00586   if (!(e.getByLabel("hcalDigis",report)))
00587     {
00588       rawOK_=false;
00589       LogWarning("HcalMonitorModule")<<" hcalDigis not available";
00590     }
00591   if (rawOK_&&!report.isValid()) {
00592     rawOK_=false;
00593   }
00594   else 
00595     {
00596       if(!fedsListed_){
00597         const std::vector<int> feds =  (*report).getFedsUnpacked();    
00598         for(unsigned int f=0; f<feds.size(); f++){
00599           meFEDS_->Fill(feds[f]);    
00600         }
00601         fedsListed_ = true;
00602       }
00603     }
00604 
00605   // try to get digis
00606   edm::Handle<HBHEDigiCollection> hbhe_digi;
00607   edm::Handle<HODigiCollection> ho_digi;
00608   edm::Handle<HFDigiCollection> hf_digi;
00609   edm::Handle<HcalTrigPrimDigiCollection> tp_digi;
00610   edm::Handle<HcalLaserDigi> laser_digi;
00611 
00612   if (!(e.getByLabel(inputLabelDigi_,hbhe_digi)))
00613     digiOK_=false;
00614   /*
00615   if (!(e.getByType(hbhe_digi)))
00616     digiOK_=false;
00617   cout <<"TEST HBHE = "<<(*hbhe_digi).size()<<endl;
00618   */
00619   if (digiOK_&&!hbhe_digi.isValid()) {
00620     digiOK_=false;
00621     LogWarning("HcalMonitorModule")<< inputLabelDigi_<<" hbhe_digi not available";
00622   }
00623 
00624   if (!(e.getByLabel(inputLabelDigi_,hf_digi)))
00625     {
00626       digiOK_=false;
00627       LogWarning("HcalMonitorModule")<< inputLabelDigi_<<" hf_digi not available";
00628     }
00629   if (digiOK_&&!hf_digi.isValid()) {
00630     digiOK_=false;
00631   }
00632 
00633   if (!(e.getByLabel(inputLabelDigi_,ho_digi)))
00634     {
00635       digiOK_=false;
00636       LogWarning("HcalMonitorModule")<< inputLabelDigi_<<" ho_digi not available";
00637     }
00638   if (digiOK_&&!ho_digi.isValid()) {
00639     digiOK_=false;
00640   }
00641   
00642   // check which Subdetectors are on by seeing which are reading out FED data
00643   // Assume subdetectors aren't present, unless we explicitly find otherwise
00644 
00645   if (digiOK_ && rawOK_)
00646     { 
00647       if ((checkHB_ && HBpresent_==0) ||
00648           (checkHE_ && HEpresent_==0) ||
00649           (checkHO_ && HOpresent_==0) ||
00650           (checkHF_ && HFpresent_==0))
00651         
00652         CheckSubdetectorStatus(*rawraw,*report,*readoutMap_,*hbhe_digi, *ho_digi, *hf_digi);
00653     }
00654   else
00655     {
00656       // Is this the behavior we want?
00657       if (debug_>1)
00658         cout <<"<HcalMonitorModule::analyze>  digiOK or rawOK error.  Assuming all subdetectors present."<<endl;
00659       HBpresent_=1;
00660       HEpresent_=1;
00661       HOpresent_=1;
00662       HFpresent_=1;
00663     }
00664 
00665   // Case where all subdetectors have no raw data -- skip event
00666   if ((checkHB_ && HBpresent_==0) &&
00667       (checkHE_ && HEpresent_==0) &&
00668       (checkHO_ && HOpresent_==0) &&
00669       (checkHF_ && HFpresent_==0))
00670     {
00671       if (debug_>1) cout <<"<HcalMonitorModule::analyze>  No HCAL raw data found for event "<<ievt_<<endl;
00672       return;
00673     }
00674 
00675   if (!(e.getByLabel(inputLabelDigi_,tp_digi)))
00676     {
00677       tpdOK_=false;
00678       LogWarning("HcalMonitorModule")<< inputLabelDigi_<<" tp_digi not available"; 
00679     }
00680 
00681   if (tpdOK_ && !tp_digi.isValid()) {
00682     tpdOK_=false;
00683   }
00684   if (!(e.getByLabel(inputLabelLaser_,laser_digi)))
00685     {laserOK_=false;}
00686   if (laserOK_&&!laser_digi.isValid()) {
00687     laserOK_=false;
00688   }
00689 
00690   // try to get rechits
00691   edm::Handle<HBHERecHitCollection> hb_hits;
00692   edm::Handle<HORecHitCollection> ho_hits;
00693   edm::Handle<HFRecHitCollection> hf_hits;
00694   edm::Handle<ZDCRecHitCollection> zdc_hits;
00695   edm::Handle<CaloTowerCollection> calotowers;
00696 
00697   if (!(e.getByLabel(inputLabelRecHitHBHE_,hb_hits)))
00698     {
00699       rechitOK_=false;
00700       //if (debug_>0)
00701         LogWarning("HcalMonitorModule")<< inputLabelRecHitHBHE_<<" not available"; 
00702     }
00703   
00704   if (rechitOK_&&!hb_hits.isValid()) {
00705     rechitOK_ = false;
00706   }
00707   if (!(e.getByLabel(inputLabelRecHitHO_,ho_hits)))
00708     {
00709       rechitOK_=false;
00710       //if (debug_>0) 
00711         LogWarning("HcalMonitorModule")<< inputLabelRecHitHO_<<" not available"; 
00712     }
00713   if (rechitOK_&&!ho_hits.isValid()) {
00714     rechitOK_ = false;
00715   }
00716   if (!(e.getByLabel(inputLabelRecHitHF_,hf_hits)))
00717     {
00718       rechitOK_=false;
00719       //if (debug_>0) 
00720         LogWarning("HcalMonitorModule")<< inputLabelRecHitHF_<<" not available"; 
00721     }
00722   if (rechitOK_&&!hf_hits.isValid()) {
00723     rechitOK_ = false;
00724   }
00725   
00726   if (!(e.getByLabel(inputLabelRecHitZDC_,zdc_hits)))
00727     {
00728       zdchitOK_=false;
00729       // ZDC Warnings should be suppressed unless debugging is on (since we don't yet normally run zdcreco)
00730       if (debug_>0) 
00731         LogWarning("HcalMonitorModule")<< inputLabelRecHitZDC_<<" not available"; 
00732     }
00733   if (zdchitOK_&&!zdc_hits.isValid()) 
00734     {
00735       zdchitOK_ = false;
00736     }
00737   
00738   // try to get calotowers 
00739   if (ctMon_!=NULL)
00740     {
00741       if (!(e.getByLabel(inputLabelCaloTower_,calotowers)))
00742         {
00743           calotowerOK_=false;
00744           if (debug_>0) LogWarning("HcalMonitorModule")<< inputLabelCaloTower_<<" not available"; 
00745         }
00746       if(calotowerOK_&&!calotowers.isValid()){
00747         calotowerOK_=false;
00748       }
00749     }
00750   else
00751     calotowerOK_=false;
00752 
00753   // Run the configured tasks, protect against missing products
00754 
00755   // Data Format monitor task
00756   
00757   if (showTiming_)
00758     {
00759       cpu_timer.reset(); cpu_timer.start();
00760     }
00761 
00762   if((dfMon_ != NULL) && (evtMask&DO_HCAL_DFMON) && rawOK_) 
00763     dfMon_->processEvent(*rawraw,*report,*readoutMap_);
00764   if (showTiming_)
00765     {
00766       cpu_timer.stop();
00767       if (dfMon_ !=NULL) cout <<"TIMER:: DATAFORMAT MONITOR ->"<<cpu_timer.cpuTime()<<endl;
00768       cpu_timer.reset(); cpu_timer.start();
00769     }
00770 
00771   if ((diTask_ != NULL) && (evtMask&DO_HCAL_DFMON) && rawOK_)
00772     diTask_->processEvent(*rawraw,*report,*readoutMap_);
00773   if (showTiming_)
00774     {
00775       cpu_timer.stop();
00776       if (diTask_ !=NULL) cout <<"TIMER:: DATA INTEGRITY TASK ->"<<cpu_timer.cpuTime()<<endl;
00777       cpu_timer.reset(); cpu_timer.start();
00778     }
00779 
00780   // Digi monitor task
00781   if((digiMon_!=NULL) && (evtMask&DO_HCAL_DIGIMON) && digiOK_) 
00782     {
00783       digiMon_->setSubDetectors(HBpresent_,HEpresent_, HOpresent_, HFpresent_);
00784       digiMon_->processEvent(*hbhe_digi,*ho_digi,*hf_digi,
00785                              *conditions_,*report);
00786     }
00787   if (showTiming_)
00788     {
00789       cpu_timer.stop();
00790       if (digiMon_ != NULL) cout <<"TIMER:: DIGI MONITOR ->"<<cpu_timer.cpuTime()<<endl;
00791       cpu_timer.reset(); cpu_timer.start();
00792     }
00793   // Pedestal monitor task
00794   if((pedMon_!=NULL) && (evtMask&DO_HCAL_PED_CALIBMON) && digiOK_) 
00795     {
00796       pedMon_->processEvent(*hbhe_digi,*ho_digi,*hf_digi,*conditions_);
00797     }
00798   if (showTiming_)
00799     {
00800       cpu_timer.stop();
00801       if (pedMon_!=NULL) cout <<"TIMER:: PEDESTAL MONITOR ->"<<cpu_timer.cpuTime()<<endl;
00802       cpu_timer.reset(); cpu_timer.start();
00803     }
00804 
00805   // LED monitor task
00806   if((ledMon_!=NULL) && (evtMask&DO_HCAL_LED_CALIBMON) && digiOK_)
00807     ledMon_->processEvent(*hbhe_digi,*ho_digi,*hf_digi,*conditions_);
00808   if (showTiming_)
00809     {
00810       cpu_timer.stop();
00811       if (ledMon_!=NULL) cout <<"TIMER:: LED MONITOR ->"<<cpu_timer.cpuTime()<<endl;
00812       cpu_timer.reset(); cpu_timer.start();
00813     }
00814 
00815   // Laser monitor task
00816   if((laserMon_!=NULL) && (evtMask&DO_HCAL_LASER_CALIBMON) && digiOK_ && laserOK_)
00817     laserMon_->processEvent(*hbhe_digi,*ho_digi,*hf_digi,*laser_digi,*conditions_);
00818   if (showTiming_)
00819     {
00820       cpu_timer.stop();
00821       if (laserMon_!=NULL) cout <<"TIMER:: LASER MONITOR ->"<<cpu_timer.cpuTime()<<endl;
00822       cpu_timer.reset(); cpu_timer.start();
00823     }
00824 
00825   // Rec Hit monitor task
00826   if((rhMon_ != NULL) && (evtMask&DO_HCAL_RECHITMON) && rechitOK_) 
00827     {
00828     rhMon_->processEvent(*hb_hits,*ho_hits,*hf_hits);
00829     // This lets us process rec hits regardless of ZDC status.
00830     // But is ZDC is okay, we'll make rec hit plots for that as well.
00831     if (zdchitOK_)
00832       {
00833         if (debug_>1) cout <<"PROCESSING ZDC!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!"<<endl;
00834         //rhMon_->processZDC(*zdc_hits);
00835       }
00836     }
00837   if (showTiming_)
00838     {
00839       cpu_timer.stop();
00840       if (rhMon_!=NULL) cout <<"TIMER:: RECHIT MONITOR ->"<<cpu_timer.cpuTime()<<endl;
00841       cpu_timer.reset(); cpu_timer.start();
00842     }
00843   
00844   // Beam Monitor task
00845   if ((beamMon_ != NULL) && (evtMask&DO_HCAL_RECHITMON) && rechitOK_)
00846     {
00847       beamMon_->processEvent(*hb_hits,*ho_hits,*hf_hits,*hf_digi);
00848     }
00849   if (showTiming_)
00850     {
00851       cpu_timer.stop();
00852       if (beamMon_!=NULL) cout <<"TIMER:: BEAM MONITOR ->"<<cpu_timer.cpuTime( \
00853 )<<endl;
00854       cpu_timer.reset(); cpu_timer.start();
00855     }
00856 
00857   // Hot Cell monitor task
00858   if((hotMon_ != NULL) && (evtMask&DO_HCAL_RECHITMON) && rechitOK_) 
00859     {
00860       hotMon_->processEvent(*hb_hits,*ho_hits,*hf_hits, 
00861                             *hbhe_digi,*ho_digi,*hf_digi,*conditions_);
00862       //hotMon_->setSubDetectors(HBpresent_,HEpresent_, HOpresent_, HFpresent_);
00863     }
00864   if (showTiming_)
00865     {
00866       cpu_timer.stop();
00867       if (hotMon_!=NULL) cout <<"TIMER:: HOTCELL MONITOR ->"<<cpu_timer.cpuTime()<<endl;
00868       cpu_timer.reset(); cpu_timer.start();
00869     }
00870   // Dead Cell monitor task -- may end up using both rec hits and digis?
00871   if((deadMon_ != NULL) && (evtMask&DO_HCAL_RECHITMON) && rechitOK_ && digiOK_) 
00872     {
00873       //deadMon_->setSubDetectors(HBpresent_,HEpresent_, HOpresent_, HFpresent_);
00874       deadMon_->processEvent(*hb_hits,*ho_hits,*hf_hits,
00875                              *hbhe_digi,*ho_digi,*hf_digi,*conditions_); 
00876     }
00877   if (showTiming_)
00878     {
00879       cpu_timer.stop();
00880       if (deadMon_!=NULL) cout <<"TIMER:: DEADCELL MONITOR ->"<<cpu_timer.cpuTime()<<endl;
00881       cpu_timer.reset(); cpu_timer.start();
00882     }
00883 
00884   // CalotowerMonitor
00885   if ((ctMon_ !=NULL) )
00886     ctMon_->processEvent(*calotowers);
00887 
00888   if (showTiming_)
00889     {
00890       cpu_timer.stop();
00891       if (ctMon_ !=NULL) cout <<"TIMER:: CALOTOWER MONITOR ->"<<cpu_timer.cpuTime()<<endl;
00892       cpu_timer.reset(); cpu_timer.start();
00893     }
00894 
00895 
00896   // Trigger Primitive monitor task -- may end up using both rec hits and digis?
00897   if((tpMon_ != NULL) && rechitOK_ && digiOK_ && tpdOK_) 
00898     tpMon_->processEvent(*hb_hits,*ho_hits,*hf_hits,
00899                          *hbhe_digi,*ho_digi,*hf_digi,*tp_digi, *readoutMap_);                       
00900   
00901 
00902   if (showTiming_)
00903     {
00904       cpu_timer.stop();
00905       if (tpMon_!=NULL) cout <<"TIMER:: TRIGGERPRIMITIVE MONITOR ->"<<cpu_timer.cpuTime()<<endl;
00906     }
00907 
00908   // Expert monitor plots
00909   if (expertMon_ != NULL) 
00910     {
00911       expertMon_->processEvent(*hb_hits,*ho_hits,*hf_hits,
00912                                *hbhe_digi,*ho_digi,*hf_digi,
00913                                *tp_digi,
00914                                *rawraw,*report,*readoutMap_);
00915     }
00916   if (showTiming_)
00917     {
00918       cpu_timer.stop();
00919       if (expertMon_!=NULL) cout <<"TIMER:: EXPERT MONITOR ->"<<cpu_timer.cpuTime()<<endl;
00920       cpu_timer.reset(); cpu_timer.start();
00921     }
00922 
00923   // Empty Event/Unsuppressed monitor plots
00924   if (eeusMon_ != NULL) 
00925     {
00926       eeusMon_->processEvent( *rawraw,*report,*readoutMap_);
00927     }
00928   if (showTiming_)
00929     {
00930       cpu_timer.stop();
00931       if (eeusMon_!=NULL) cout <<"TIMER:: EE/US MONITOR ->"<<cpu_timer.cpuTime()<<endl;
00932       cpu_timer.reset(); cpu_timer.start();
00933     }
00934 
00935 
00936   if(debug_>0 && ievt_%1000 == 0)
00937     cout << "HcalMonitorModule: processed " << ievt_ << " events" << endl;
00938 
00939   if(debug_>1)
00940     {
00941       cout << "HcalMonitorModule: processed " << ievt_ << " events" << endl;
00942       cout << "    RAW Data   ==> " << rawOK_<< endl;
00943       cout << "    Digis      ==> " << digiOK_<< endl;
00944       cout << "    RecHits    ==> " << rechitOK_<< endl;
00945       cout << "    TrigRec    ==> " << trigOK_<< endl;
00946       cout << "    TPdigis    ==> " << tpdOK_<< endl;    
00947       cout << "    CaloTower  ==> " << calotowerOK_ <<endl;
00948       cout << "    LaserDigis ==> " << laserOK_ << endl;
00949     }
00950   
00951   return;
00952 }
00953 
00954 //--------------------------------------------------------
00955 bool HcalMonitorModule::prescale()
00956 {
00959   if (debug_>0) cout <<"HcalMonitorModule::prescale"<<endl;
00960   
00961   gettimeofday(&psTime_.updateTV,NULL);
00962   double time = (psTime_.updateTV.tv_sec*1000.0+psTime_.updateTV.tv_usec/1000.0);
00963   time/= (1000.0); 
00964   psTime_.elapsedTime = time - psTime_.updateTime;
00965   psTime_.updateTime = time;
00966   //First determine if we care...
00967   bool evtPS =    prescaleEvt_>0;
00968   bool lsPS =     prescaleLS_>0;
00969   bool timePS =   prescaleTime_>0;
00970   bool updatePS = prescaleUpdate_>0;
00971 
00972   // If no prescales are set, keep the event
00973   if(!evtPS && !lsPS && !timePS && !updatePS)
00974     {
00975       return false;
00976     }
00977   //check each instance
00978   if(lsPS && (ilumisec_%prescaleLS_)!=0) lsPS = false; //LS veto
00979   //if(evtPS && (ievent_%prescaleEvt_)!=0) evtPS = false; //evt # veto
00980   // we can't just call (ievent_%prescaleEvt_) because ievent values not consecutive
00981   if (evtPS && (ievt_pre_%prescaleEvt_)!=0) evtPS = false;
00982   if(timePS)
00983     {
00984       double elapsed = (psTime_.updateTime - psTime_.vetoTime)/60.0;
00985       if(elapsed<prescaleTime_){
00986         timePS = false;  //timestamp veto
00987         psTime_.vetoTime = psTime_.updateTime;
00988       }
00989     } //if (timePS)
00990 
00991   //  if(prescaleUpdate_>0 && (nupdates_%prescaleUpdate_)==0) updatePS=false; ///need to define what "updates" means
00992   
00993   if (debug_>1) 
00994     {
00995       cout<<"HcalMonitorModule::prescale  evt: "<<ievent_<<"/"<<evtPS<<", ";
00996       cout <<"ls: "<<ilumisec_<<"/"<<lsPS<<",";
00997       cout <<"time: "<<psTime_.updateTime - psTime_.vetoTime<<"/"<<timePS<<endl;
00998     }  
00999   // if any criteria wants to keep the event, do so
01000   if(evtPS || lsPS || timePS) return false; //FIXME updatePS left out for now
01001   return true;
01002 } // HcalMonitorModule::prescale(...)
01003 
01004 
01005 void HcalMonitorModule::CheckSubdetectorStatus(const FEDRawDataCollection& rawraw, 
01006                                                const HcalUnpackerReport& report, 
01007                                                const HcalElectronicsMap& emap,
01008                                                const HBHEDigiCollection& hbhedigi,
01009                                                const HODigiCollection& hodigi,
01010                                                const HFDigiCollection& hfdigi
01011                                                //const ZDCDigiCollection& zdcdigi,
01012 
01013                                                )
01014 {
01015   vector<int> fedUnpackList;
01016   for (int i=FEDNumbering::getHcalFEDIds().first; 
01017        i<=FEDNumbering::getHcalFEDIds().second; 
01018        i++) 
01019     {
01020       fedUnpackList.push_back(i);
01021     }
01022   
01023   for (vector<int>::const_iterator i=fedUnpackList.begin();
01024        i!=fedUnpackList.end(); 
01025        ++i) 
01026     {
01027       const FEDRawData& fed = rawraw.FEDData(*i);
01028       if (fed.size()<12) continue; // Was 16. How do such tiny events even get here?
01029       
01030       // get the DCC header
01031       const HcalDCCHeader* dccHeader=(const HcalDCCHeader*)(fed.data());
01032       if (!dccHeader) return;
01033       int dccid=dccHeader->getSourceId();
01034       // check for HF
01035       if (dccid>717 && dccid<724)
01036         {
01037           if (HFpresent_==0 && hfdigi.size()>0)
01038             {
01039               HFpresent_ = 1;
01040               meHF_->Fill(HFpresent_);
01041             }
01042           continue;
01043         }
01044 
01045       // check for HO
01046       if (dccid>723)
01047         {
01048           if (HOpresent_==0 && hodigi.size()>0)
01049             {
01050               HOpresent_ = 1;
01051               meHO_->Fill(HOpresent_);
01052             }
01053           continue;
01054         }
01055       
01056       // Looking at HB and HE is more complicated, since they're combined into HBHE
01057       // walk through the HTR data...
01058       HcalHTRData htr;  
01059       for (int spigot=0; spigot<HcalDCCHeader::SPIGOT_COUNT; spigot++) {    
01060         if (!dccHeader->getSpigotPresent(spigot)) continue;
01061         
01062         // Load the given decoder with the pointer and length from this spigot.
01063         dccHeader->getSpigotData(spigot,htr, fed.size()); 
01064         
01065         // check min length, correct wordcount, empty event, or total length if histo event.
01066         if (!htr.check()) continue;
01067         if (htr.isHistogramEvent()) continue;
01068         
01069         int firstFED =  FEDNumbering::getHcalFEDIds().first; 
01070         
01071         // Tease out HB and HE, which share HTRs in HBHE
01072         for(int fchan=0; fchan<3; ++fchan) //0,1,2 are valid
01073           {
01074             for(int fib=1; fib<9; ++fib) //1...8 are valid
01075               {
01076                 HcalElectronicsId eid(fchan,fib,spigot,dccid-firstFED);
01077                 eid.setHTR(htr.readoutVMECrateId(),
01078                            htr.htrSlot(),htr.htrTopBottom());
01079                 DetId did=emap.lookup(eid);
01080                 if (!did.null()) 
01081                   {
01082                     
01083                     switch (((HcalSubdetector)did.subdetId()))
01084                       {
01085                       case (HcalBarrel): 
01086                         {
01087                           if (HBpresent_==0)
01088                             {
01089                               HBpresent_ = 1;
01090                               meHB_->Fill(HBpresent_);
01091                             }
01092                         } break; // case (HcalBarrel)
01093                       case (HcalEndcap): 
01094                         {
01095                           if (HEpresent_==0)
01096                             {
01097                               HEpresent_ = 1;
01098                               meHE_->Fill(HEpresent_);
01099                             }
01100                         } break; // case (HcalEndcap)
01101                       case (HcalOuter): 
01102                         { // shouldn't reach these last two cases
01103                           if (HOpresent_==0)
01104                             {
01105                               {
01106                                 HOpresent_ = 1;
01107                                 meHO_->Fill(HOpresent_);
01108                                 return;
01109                               }
01110                             } 
01111                         } break; // case (HcalOuter)
01112                       case (HcalForward): 
01113                         {
01114                           if (HFpresent_==0)
01115                             {
01116                               meHF_->Fill(HFpresent_);
01117                               HFpresent_ = 1;
01118                             }
01119                         } break; //case (HcalForward)
01120                       default: break;
01121                       } // switch ((HcalSubdetector...)
01122                   } // if (!did.null())
01123               } // for (int fib=0;...)
01124           } // for (int fchan = 0;...)
01125         
01126       } // for (int spigot=0;...)
01127     } //  for (vector<int>::const_iterator i=fedUnpackList.begin();
01128   return;
01129 } // void HcalMonitorModule::CheckSubdetectorStatus(...)
01130 
01131 #include "FWCore/Framework/interface/MakerMacros.h"
01132 #include <DQM/HcalMonitorModule/src/HcalMonitorModule.h>
01133 #include "DQMServices/Core/interface/DQMStore.h"
01134 
01135 DEFINE_FWK_MODULE(HcalMonitorModule);

Generated on Tue Jun 9 17:32:57 2009 for CMSSW by  doxygen 1.5.4