CMS 3D CMS Logo

RPCMonitorDigi Class Reference

#include <DQM/RPCMonitorDigi/interface/RPCMonitorDigi.h>

Inheritance diagram for RPCMonitorDigi:

edm::EDAnalyzer

List of all members.

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
virtual void beginJob (edm::EventSetup const &)
void beginRun (const edm::Run &r, const edm::EventSetup &c)
std::map< std::string,
MonitorElement * > 
bookDetUnitME (RPCDetId &, const edm::EventSetup &)
 Booking of MonitoringElemnt for one RPCDetId (= roll).
std::map< std::string,
MonitorElement * > 
bookRegionRing (int region, int ring)
 Booking of MonitoringElemnt at Wheel/Disk level.
virtual void endJob (void)
 RPCMonitorDigi (const edm::ParameterSet &)
 ~RPCMonitorDigi ()

Private Member Functions

int stripsInRoll (RPCDetId &, const edm::EventSetup &)

Private Attributes

MonitorElementBarrelOccupancy
MonitorElementClusterSize_for_Barrel
MonitorElementClusterSize_for_BarrelandEndcaps
MonitorElementClusterSize_for_EndcapNegative
MonitorElementClusterSize_for_EndcapPositive
int counter
DQMStoredbe
 DQM store.
bool dqmexpert
bool dqmshifter
bool dqmsuperexpert
MonitorElementEndcapNegativeOccupancy
MonitorElementEndcapPositiveOccupancy
std::map< uint32_t, boolfoundHitsInChamber
std::string GlobalHistogramsFolder
std::map< uint32_t, std::map
< std::string, MonitorElement * > > 
meCollection
bool mergeRuns_
std::map< std::pair< int, int >,
std::map< std::string,
MonitorElement * > > 
meWheelDisk
std::string nameInLog
MonitorElementNumberOfClusters_for_Barrel
MonitorElementNumberOfClusters_for_EndcapNegative
MonitorElementNumberOfClusters_for_EndcapPositive
MonitorElementNumberOfDigis_for_Barrel
MonitorElementNumberOfDigis_for_EndcapNegative
MonitorElementNumberOfDigis_for_EndcapPositive
std::string RootFileName
MonitorElementSameBxDigisMeBarrel_
MonitorElementSameBxDigisMeEndcapNegative_
MonitorElementSameBxDigisMeEndcapPositive_
bool saveRootFile
int saveRootFileEventsInterval


Detailed Description

Definition at line 32 of file RPCMonitorDigi.h.


Constructor & Destructor Documentation

RPCMonitorDigi::RPCMonitorDigi ( const edm::ParameterSet pset  )  [explicit]

Definition at line 29 of file RPCMonitorDigi.cc.

References dqmexpert, dqmshifter, dqmsuperexpert, foundHitsInChamber, edm::ParameterSet::getUntrackedParameter(), mergeRuns_, nameInLog, RootFileName, saveRootFile, and saveRootFileEventsInterval.

00029                                                         :counter(0){
00030   foundHitsInChamber.clear();
00031   nameInLog = pset.getUntrackedParameter<string>("moduleLogName", "RPC_DQM");
00032 
00033   saveRootFile  = pset.getUntrackedParameter<bool>("DigiDQMSaveRootFile", false); 
00034   mergeRuns_  = pset.getUntrackedParameter<bool>("MergeDifferentRuns", false); 
00035   saveRootFileEventsInterval  = pset.getUntrackedParameter<int>("DigiEventsInterval", 10000);
00036   RootFileName  = pset.getUntrackedParameter<string>("RootFileNameDigi", "RPCMonitor.root"); 
00037 
00038   dqmshifter = pset.getUntrackedParameter<bool>("dqmshifter", false);
00039   dqmexpert = pset.getUntrackedParameter<bool>("dqmexpert", false);
00040   dqmsuperexpert = pset.getUntrackedParameter<bool>("dqmsuperexpert", false);
00041 }

RPCMonitorDigi::~RPCMonitorDigi (  ) 

Definition at line 43 of file RPCMonitorDigi.cc.

00043 {}


Member Function Documentation

void RPCMonitorDigi::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
) [virtual]

RPC Geometry

DIGI

RecHits

end loop on RPCRecHits for given roll

end loop on RPC Digi Collection

must be fixed!!!!!!!!!!!!!!!!!!!!!!!!!!!!!

Implements edm::EDAnalyzer.

Definition at line 155 of file RPCMonitorDigi.cc.

References BarrelOccupancy, bookDetUnitME(), bookRegionRing(), ClusterSize_for_Barrel, ClusterSize_for_BarrelandEndcaps, ClusterSize_for_EndcapNegative, ClusterSize_for_EndcapPositive, counter, detId, rpcdqm::utils::detId2RollNr(), dqmexpert, dqmsuperexpert, EndcapNegativeOccupancy, EndcapPositiveOccupancy, error, MonitorElement::Fill(), find(), foundHitsInChamber, edm::EventSetup::get(), edm::Event::getByType(), index, it, k, meCollection, meWheelDisk, VarParsing::mult, RPCGeomServ::name(), nameInLog, NumberOfClusters_for_Barrel, NumberOfClusters_for_EndcapNegative, NumberOfClusters_for_EndcapPositive, NumberOfDigis_for_Barrel, NumberOfDigis_for_EndcapNegative, NumberOfDigis_for_EndcapPositive, RPCDetId::region(), RPCDetId::ring(), SameBxDigisMeBarrel_, RPCDetId::sector(), RPCDetId::station(), strip(), GeomDet::surface(), Surface::toGlobal(), PV3DBase< T, PVType, FrameType >::x(), x, LocalError::xx(), and PV3DBase< T, PVType, FrameType >::y().

00155                                                                          {
00156   counter++;
00157   LogInfo (nameInLog) <<"[RPCMonitorDigi]: Beginning analyzing event " << counter;
00158   
00159   map<uint32_t, bool >::iterator mapItrReset;
00160   for (mapItrReset = foundHitsInChamber.begin(); mapItrReset != foundHitsInChamber.end(); ++ mapItrReset) {
00161     mapItrReset->second=false;
00162   }
00163   
00165   ESHandle<RPCGeometry> rpcGeo;
00166   iSetup.get<MuonGeometryRecord>().get(rpcGeo);
00167   
00169   Handle<RPCDigiCollection> rpcdigis;
00170   iEvent.getByType(rpcdigis);
00171 
00173   Handle<RPCRecHitCollection> rpcHits;
00174   iEvent.getByType(rpcHits);
00175 
00176   map<int,int> bxMap;
00177   int totalNumberDigi[3][1];
00178 
00179   for (int k=1; k<=3;k++){
00180     totalNumberDigi[k][1]=0;
00181   }
00182 
00183   RPCDigiCollection::DigiRangeIterator collectionItr;
00184   //Loop on digi collection
00185   for(collectionItr=rpcdigis->begin(); collectionItr!=rpcdigis->end(); ++collectionItr){
00186   
00187     RPCDetId detId=(*collectionItr ).first; 
00188     uint32_t id=detId(); 
00189 
00190     const GeomDet* gdet=rpcGeo->idToDet(detId);
00191     const BoundPlane & surface = gdet->surface();
00192     
00193     //get roll name
00194     RPCGeomServ RPCname(detId);
00195     string nameRoll = RPCname.name();
00196  
00197     stringstream os;
00198 
00199     rpcdqm::utils prova;
00200     int nr = prova.detId2RollNr(detId);
00201 
00202     std::map<uint32_t, std::map<std::string,MonitorElement*> >::iterator meItr = meCollection.find(id);
00203   
00204     if (meItr == meCollection.end() || (meCollection.size()==0)) {
00205       meCollection[id]=bookDetUnitME(detId,iSetup );
00206     }
00207   
00208     std::map<std::string, MonitorElement*> meMap=meCollection[id];
00209   
00210     int region=detId.region();
00211     int ring;
00212     string regionName;
00213     string ringType;
00214     if(detId.region() == 0) {
00215       regionName="Barrel";  
00216       ringType = "Wheel";  
00217       ring = detId.ring();
00218     }else if (detId.region() == -1){
00219       regionName="Encap-";
00220       ringType =  "Disk";
00221       ring = detId.region()*detId.station();
00222     }else{
00223       regionName="Encap+";
00224       ringType =  "Disk";
00225       ring = detId.station();
00226     }
00227     std::pair<int,int> regionRing(region,ring);
00228     std::map<std::pair<int,int>, std::map<std::string,MonitorElement*> >::iterator meRingItr = meWheelDisk.find(regionRing);
00229     if (meRingItr == meWheelDisk.end() || (meWheelDisk.size()==0)) {
00230       meWheelDisk[regionRing]=bookRegionRing(region,ring);
00231     }
00232  
00233     map<std::string, MonitorElement*> meRingMap=meWheelDisk[regionRing];
00234  
00235     vector<int> strips;
00236     vector<int> bxs;
00237     strips.clear(); 
00238     bxs.clear(); 
00239 
00240     //get the RecHits associated to the roll
00241     typedef pair<RPCRecHitCollection::const_iterator, RPCRecHitCollection::const_iterator> rangeRecHits;
00242     rangeRecHits recHitCollection =  rpcHits->get(detId);
00243  
00244     int numberOfDigi= 0;
00245  
00246     //remove duplicated digis
00247     RPCDigiCollection::const_iterator digiItr; 
00248     //loop on digis of given roll
00249      for (digiItr =(*collectionItr ).second.first;digiItr != (*collectionItr ).second.second; ++digiItr){
00250       int strip= (*digiItr).strip();
00251       int bx=(*digiItr).bx();
00252       // pair<int, int> stripBx(strip,bx);
00253      
00254       //skip 
00255       vector<int>::const_iterator itrStrips = find(strips.begin(),strips.end(),strip);
00256       if(itrStrips!=strips.end() && strips.size()!=0) continue;
00257       ++numberOfDigi;
00258       strips.push_back(strip);
00259       
00260       //get bx number for this digi
00261       vector<int>::const_iterator existingBX = find(bxs.begin(),bxs.end(),bx);
00262       if(existingBX==bxs.end())bxs.push_back(bx);
00263    
00264       //adding new histo C.Carrillo & A. Cimmino
00265       map<int,int>::const_iterator bxItr = bxMap.find((*digiItr).bx());
00266       if (bxItr == bxMap.end()|| bxMap.size()==0 )bxMap[(*digiItr).bx()]=1;
00267       else bxMap[(*digiItr).bx()]++;
00268    
00269       //sector based histograms for dqm shifter
00270       os.str("");
00271       os<<"1DOccupancy_"<<ringType<<"_"<<ring;
00272       string meId = os.str();
00273       meRingMap[meId]->Fill(detId.sector());
00274       os.str("");
00275       os<<"Sec"<<detId.sector();
00276       meRingMap[meId] -> setBinLabel(detId.sector(), os.str(), 1);
00277    
00278       os.str("");
00279       os<<"BxDistribution_"<<ringType<<"_"<<ring<<"_Sector_"<<detId.sector();
00280       meMap[os.str()]->Fill(bx);
00281    
00282       os.str("");
00283       os<<"BxDistribution_"<<ringType<<"_"<<ring;
00284       meRingMap[os.str()]->Fill(bx);
00285    
00286       if(detId.region()==0)
00287         BarrelOccupancy -> Fill(detId.sector(), ring);
00288       else if(detId.region()==1)
00289         EndcapPositiveOccupancy -> Fill(detId.sector(), ring);
00290       else if(detId.region()==-1)
00291         EndcapNegativeOccupancy -> Fill(detId.sector(), ring);
00292 
00293       os.str("");
00294       os<<"Occupancy_"<<ringType<<"_"<<ring<<"_Sector_"<<detId.sector();
00295       meMap[os.str()]->Fill(strip, nr);
00296           
00297       string Yaxis=nameRoll;
00298       if (detId.region()==0){
00299         Yaxis.erase (1,1);
00300         Yaxis.erase(0,3);
00301         Yaxis.replace(Yaxis.find("S"),4,"");
00302         Yaxis.erase(Yaxis.find("_")+2,8);
00303       }else{
00304         Yaxis.erase(0,8);
00305       }
00306       meMap[os.str()]->setBinLabel(nr, Yaxis, 2);
00307       
00308       os.str("");
00309       os<<"Occupancy_"<<nameRoll;
00310       meMap[os.str()]->Fill(strip);
00311   
00312       if(dqmexpert){    
00313         os.str("");
00314         os<<"BXN_"<<nameRoll;
00315         meMap[os.str()]->Fill(bx);
00316         }
00317   
00318       if (dqmsuperexpert) {     
00319         os.str("");
00320         os<<"BXN_vs_strip_"<<nameRoll;
00321         meMap[os.str()]->Fill(strip,bx);
00322       }
00323     }  //end loop of digis of given roll
00324   
00325     if (dqmexpert){
00326 
00327       for(unsigned int stripIter=0;stripIter<strips.size(); ++stripIter){       
00328         if(stripIter< (strips.size()-1) && strips[stripIter+1]==strips[stripIter]+1) {
00329           os.str("");
00330           os<<"CrossTalkHigh_"<<nameRoll;
00331           meMap[os.str()]->Fill(strips[stripIter]);     
00332         }
00333         if(stripIter >0 && strips[stripIter-1]==strips[stripIter]-1) {
00334           os.str("");
00335           os<<"CrossTalkLow_"<<nameRoll;
00336               meMap[os.str()]->Fill(strips[stripIter]); 
00337         }
00338       }
00339       
00340       os.str("");
00341       os<<"NumberOfDigi_"<<nameRoll;
00342       meMap[os.str()]->Fill(numberOfDigi);
00343       
00344       os.str("");
00345       os<<"BXWithData_"<<nameRoll;
00346       meMap[os.str()]->Fill(bxs.size());
00347     }
00348  
00349     os.str("");
00350     os<<"BXWithData_"<<ringType<<"_"<<ring<<"_Sector_"<<detId.sector();
00351     meMap[os.str()]->Fill(bxs.size());
00352  
00353     os.str("");
00354     os<<"NumberOfDigi_"<<ringType<<"_"<<ring<<"_Sector_"<<detId.sector();
00355     meMap[os.str()]->Fill(numberOfDigi);
00356  
00357     totalNumberDigi[detId.region()+2][1]+= numberOfDigi;        
00358  
00359     // Fill RecHit MEs   
00360     if(recHitCollection.first==recHitCollection.second ){   
00361  
00362       if(dqmsuperexpert) {
00363         os.str("");
00364         os<<"MissingHits_"<<nameRoll;
00365         meMap[os.str()]->Fill((int)(counter), 1.0);
00366       }
00367     }else{     
00368       foundHitsInChamber[id]=true;
00369  
00370       RPCRecHitCollection::const_iterator it;
00371       int numberOfHits=0;    
00372       int numbOfClusters=0;
00373       //loop RPCRecHits for given roll
00374       for (it = recHitCollection.first; it != recHitCollection.second ; it++) {
00375         numbOfClusters++; 
00376         RPCDetId detIdRecHits=it->rpcId();
00377         LocalError error=it->localPositionError();//plot of errors/roll => should be gaussian   
00378         LocalPoint point=it->localPosition();     //plot of coordinates/roll =>should be flat
00379         GlobalPoint globalHitPoint=surface.toGlobal(point); 
00380  
00381         os.str("");
00382         os<<"OccupancyXY_"<<ringType<<"_"<<ring;
00383         meRingMap[os.str()]->Fill(globalHitPoint.x(),globalHitPoint.y());
00384 
00385         int mult=it->clusterSize();               //cluster size plot => should be within 1-3   
00386         int firstStrip=it->firstClusterStrip();    //plot first Strip => should be flat
00387         float xposition=point.x();
00388 
00389         ClusterSize_for_BarrelandEndcaps -> Fill(mult);
00390 
00391         if(detId.region() ==  0) {
00392           ClusterSize_for_Barrel -> Fill(mult);
00393         } else if (detId.region() ==  -1) {
00394           if(mult<=10) ClusterSize_for_EndcapNegative -> Fill(mult);
00395           else ClusterSize_for_EndcapNegative -> Fill(11);         
00396         } else if (detId.region() ==  1) {
00397           if(mult<=10) ClusterSize_for_EndcapPositive -> Fill(mult);
00398           else ClusterSize_for_EndcapPositive -> Fill(11);
00399         } 
00400 
00401         //Cluster Size by Wheels and sector
00402         os.str("");
00403         os<<"ClusterSize_"<<ringType<<"_"<<ring;
00404         meRingMap[os.str()] -> Fill(mult); 
00405 
00406         os.str("");
00407         os<<"ClusterSize_"<<ringType<<"_"<<ring<<"_Sector_"<<detId.sector();
00408         meMap[os.str()] -> Fill(mult);
00409 
00410         if (dqmsuperexpert) {
00411           int centralStrip=firstStrip;
00412           if(mult%2) {
00413             centralStrip+= mult/2;
00414           }else{        
00415             float x = gRandom->Uniform(2);
00416             centralStrip+=(x<1)? (mult/2)-1 : (mult/2);
00417           }
00418           os.str("");
00419           os<<"ClusterSize_vs_CentralStrip_"<<nameRoll;
00420           meMap[os.str()]->Fill(centralStrip,mult);
00421         
00422           for(int index=0; index<mult; ++index){
00423             os.str("");
00424             os<<"ClusterSize_vs_Strip_"<<nameRoll;
00425             meMap[os.str()]->Fill(firstStrip+index,mult);
00426           }
00427 
00428           os.str("");
00429           os<<"ClusterSize_vs_LowerSrip_"<<nameRoll;
00430           meMap[os.str()]->Fill(firstStrip,mult);
00431         
00432           os.str("");
00433           os<<"ClusterSize_vs_HigherStrip_"<<nameRoll;
00434           meMap[os.str()]->Fill(firstStrip+mult-1,mult);
00435         
00436           os.str("");
00437           os<<"RecHitX_vs_dx_"<<nameRoll;
00438           meMap[os.str()]->Fill(xposition,error.xx());
00439         }
00440 
00441         if(dqmexpert) {
00442           os.str("");
00443           os<<"ClusterSize_"<<nameRoll;
00444           meMap[os.str()]->Fill(mult);
00445           
00446           os.str("");
00447           os<<"RecHitXPosition_"<<nameRoll;
00448           meMap[os.str()]->Fill(xposition);
00449           
00450           os.str("");
00451           os<<"RecHitDX_"<<nameRoll;
00452           meMap[os.str()]->Fill(error.xx());       
00453         }
00454         numberOfHits++;
00455       }
00456       
00457 
00458       if(dqmexpert) {    
00459         os.str("");
00460         os<<"NumberOfClusters_"<<nameRoll;
00461         meMap[os.str()]->Fill(numbOfClusters);
00462         
00463         if(numberOfHits>5) numberOfHits=16;
00464         os.str("");
00465         os<<"RecHitCounter_"<<nameRoll;
00466         meMap[os.str()]->Fill(numberOfHits);
00467       }
00468       
00469       if(detId.region()==0)
00470         NumberOfClusters_for_Barrel -> Fill(numbOfClusters);
00471       else if (detId.region()==1)
00472         NumberOfClusters_for_EndcapPositive -> Fill(numbOfClusters);
00473       else if(detId.region()==-1)
00474         NumberOfClusters_for_EndcapNegative -> Fill(numbOfClusters);
00475       os.str("");
00476       os<<"NumberOfClusters_"<<ringType<<"_"<<ring<<"_Sector_"<<detId.sector();
00477       meMap[os.str()]->Fill(numbOfClusters);
00478     }
00479   }
00480  
00481   //fill global histo at the end of the event
00482       NumberOfDigis_for_Barrel ->Fill(totalNumberDigi[2][1]);
00483       NumberOfDigis_for_EndcapPositive ->Fill(totalNumberDigi[3][1]);
00484       NumberOfDigis_for_EndcapNegative ->Fill(totalNumberDigi[1][1]);
00485 
00486   //adding new histo C.Carrillo & A. Cimmino
00487   for (map<int, int>::const_iterator myItr= bxMap.begin(); 
00488        myItr!=bxMap.end(); myItr++){
00489     SameBxDigisMeBarrel_ ->Fill((*myItr).second);
00490   } 
00491 }

void RPCMonitorDigi::beginJob ( edm::EventSetup const &   )  [virtual]

get hold of back-end interface

Reimplemented from edm::EDAnalyzer.

Definition at line 45 of file RPCMonitorDigi.cc.

References BarrelOccupancy, DQMStore::book1D(), ClusterSize_for_Barrel, ClusterSize_for_BarrelandEndcaps, ClusterSize_for_EndcapNegative, ClusterSize_for_EndcapPositive, dbe, EndcapNegativeOccupancy, EndcapPositiveOccupancy, GlobalHistogramsFolder, i, nameInLog, NumberOfClusters_for_Barrel, NumberOfClusters_for_EndcapNegative, NumberOfClusters_for_EndcapPositive, NumberOfDigis_for_Barrel, NumberOfDigis_for_EndcapNegative, NumberOfDigis_for_EndcapPositive, SameBxDigisMeBarrel_, and DQMStore::setCurrentFolder().

00045                                               {
00046   LogInfo (nameInLog) <<"[RPCMonitorDigi]: Begin job" ;
00047   
00049   dbe = Service<DQMStore>().operator->();
00050 
00051   GlobalHistogramsFolder="RPC/RecHits/SummaryHistograms";
00052   dbe->setCurrentFolder(GlobalHistogramsFolder);  
00053 
00054   ClusterSize_for_Barrel = dbe->book1D("ClusterSize_for_Barrel", "ClusterSize for Barrel", 20, 0.5, 20.5);
00055   ClusterSize_for_EndcapPositive = dbe->book1D("ClusterSize_for_EndcapPositive", "ClusterSize for PositiveEndcap",  20, 0.5, 20.5);
00056   ClusterSize_for_EndcapNegative = dbe->book1D("ClusterSize_for_EndcapNegative", "ClusterSize for NegativeEndcap", 20, 0.5, 20.5);
00057 
00058   ClusterSize_for_BarrelandEndcaps = dbe->book1D("ClusterSize_for_BarrelandEndcap", "ClusterSize for Barrel&Endcaps", 20, 0.5, 20.5);
00059 
00060   NumberOfDigis_for_Barrel = dbe -> book1D("NumberOfDigi_for_Barrel", "NumberOfDigis for Barrel", 200, 0.5, 200.5);
00061   NumberOfDigis_for_EndcapPositive = dbe -> book1D("NumberOfDigi_for_EndcapPositive", "NumberOfDigis for Endcap Positive", 200, 0.5, 200.5);
00062   NumberOfDigis_for_EndcapNegative = dbe -> book1D("NumberOfDigi_for_EndcapNegative", "NumberOfDigis for Endcap Negative", 200, 0.5, 200.5);
00063 
00064   NumberOfClusters_for_Barrel = dbe -> book1D("NumberOfClusters_for_Barrel", "NumberOfClusters for Barrel", 20, 0.5, 20.5);
00065   NumberOfClusters_for_EndcapPositive = dbe -> book1D("NumberOfClusters_for_EndcapPositive", "NumberOfClusters for Endcap Positive", 20, 0.5, 20.5);
00066   NumberOfClusters_for_EndcapNegative = dbe -> book1D("NumberOfClusters_for_EndcapNegative", "NumberOfClusters for Endcap Negative", 20, 0.5, 20.5);
00067 
00068   SameBxDigisMeBarrel_ = dbe->book1D("SameBXDigis_Barrel", "Digis with same bx", 20, 0.5, 20.5);  
00069  //  SameBxDigisMeEndcapPositive_ = dbe->book1D("SameBXDigis_EndcapPositive", "Digis with same bx", 20, 0.5, 20.5);  
00070  //   SameBxDigisMeEndcapNegative_ = dbe->book1D("SameBXDigis_EndcapNegative", "Digis with same bx", 20, 0.5, 20.5);  
00071 
00072   BarrelOccupancy = dbe -> book2D("Occupancy_for_Barrel", "Barrel Occupancy Wheel vs Sector", 12, 0.5, 12.5, 5, -2.5, 2.5);
00073   EndcapPositiveOccupancy = dbe -> book2D("Occupancy_for_EndcapPositive", "Endcap Positive Occupancy Disk vs Sector", 6, 0.5, 6.5, 8, -4, 4);
00074   EndcapNegativeOccupancy = dbe -> book2D("Occupancy_for_EndcapNegative", "Endcap Negative Occupancy Disk vs Sector", 6, 0.5, 6.5, 8, -4, 4);
00075  
00076   stringstream binLabel;
00077   for (int i = 1; i<13; i++){
00078     binLabel.str("");
00079     binLabel<<"Sec"<<i;
00080     BarrelOccupancy -> setBinLabel(i, binLabel.str(), 1);
00081     if(i<6){
00082       binLabel.str("");
00083       binLabel<<"Wheel"<<i-3;
00084       BarrelOccupancy -> setBinLabel(i, binLabel.str(), 2);
00085     }    
00086     if(i<7) {
00087       binLabel.str("");
00088       binLabel<<"Sec"<<i;
00089       EndcapPositiveOccupancy -> setBinLabel(i, binLabel.str(), 1);
00090       EndcapNegativeOccupancy -> setBinLabel(i, binLabel.str(), 1);
00091     }
00092       if(i<9){
00093       binLabel.str("");
00094       if (i<5) binLabel<<"Disk"<<(i-1)-4; else  binLabel<<"Disk"<<i-4;
00095       EndcapPositiveOccupancy -> setBinLabel(i, binLabel.str(), 2);
00096       EndcapNegativeOccupancy -> setBinLabel(i, binLabel.str(), 2);
00097     }
00098   }
00099 }

void RPCMonitorDigi::beginRun ( const edm::Run r,
const edm::EventSetup c 
) [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 101 of file RPCMonitorDigi.cc.

References BarrelOccupancy, ClusterSize_for_Barrel, ClusterSize_for_BarrelandEndcaps, ClusterSize_for_EndcapNegative, ClusterSize_for_EndcapPositive, EndcapNegativeOccupancy, EndcapPositiveOccupancy, meCollection, mergeRuns_, meWheelDisk, nameInLog, NumberOfClusters_for_Barrel, NumberOfClusters_for_EndcapNegative, NumberOfClusters_for_EndcapPositive, NumberOfDigis_for_Barrel, NumberOfDigis_for_EndcapNegative, NumberOfDigis_for_EndcapPositive, MonitorElement::Reset(), and SameBxDigisMeBarrel_.

00101                                                               {
00102   LogInfo (nameInLog) <<"Begin Run " ;
00103   //if mergeRuns_ skip reset
00104   //if merge remember to safe at job end and not at run end
00105   if (mergeRuns_) return;
00106 
00107   //MEs are reset at every new run. They are saved at the end of each run
00108   //Reset all Histograms
00109   for (map<uint32_t, map<string,MonitorElement*> >::const_iterator meItr = meCollection.begin();
00110        meItr!= meCollection.end(); meItr++){
00111     for (map<string,MonitorElement*>::const_iterator Itr = (*meItr).second.begin();
00112          Itr!= (*meItr).second.end(); Itr++){
00113       (*Itr).second->Reset();
00114     }
00115   }
00116 
00117   for (map<pair<int,int>, map<string, MonitorElement*> >::const_iterator meItr =  meWheelDisk.begin();
00118        meItr!=  meWheelDisk.end(); meItr++){
00119     for (map<string,MonitorElement*>::const_iterator Itr = (*meItr).second.begin();
00120          Itr!= (*meItr).second.end(); Itr++){
00121       (*Itr).second->Reset();
00122     }
00123   }
00124   
00125   //Reset All Global histos !!!!!!!!!!!!!!!!!!!!!!
00126   ClusterSize_for_Barrel->Reset();
00127   ClusterSize_for_EndcapPositive ->Reset();
00128   ClusterSize_for_EndcapNegative->Reset();
00129   ClusterSize_for_BarrelandEndcaps->Reset(); 
00130   
00131   NumberOfDigis_for_Barrel ->Reset();
00132   NumberOfDigis_for_EndcapPositive ->Reset();
00133   NumberOfDigis_for_EndcapNegative->Reset();
00134 
00135   NumberOfClusters_for_Barrel ->Reset();
00136   NumberOfClusters_for_EndcapPositive->Reset(); 
00137   NumberOfClusters_for_EndcapNegative ->Reset();
00138 
00139 
00140   SameBxDigisMeBarrel_->Reset();
00141  //  SameBxDigisMeEndcapPositive_->Reset();
00142 //   SameBxDigisMeEndcapNegative_ ->Reset();
00143   
00144   BarrelOccupancy ->Reset();
00145   EndcapPositiveOccupancy ->Reset();
00146   EndcapNegativeOccupancy->Reset();
00147 }

map< string, MonitorElement * > RPCMonitorDigi::bookDetUnitME ( RPCDetId detId,
const edm::EventSetup iSetup 
)

Booking of MonitoringElemnt for one RPCDetId (= roll).

Name components common to current RPCDetId

RPCRecHits

RPCRecHits

Definition at line 10 of file RPCBookDetUnitME.cc.

References DQMStore::book1D(), DQMStore::book2D(), dbe, dqmexpert, dqmsuperexpert, RPCBookFolderStructure::folderStructure(), DQMStore::get(), RPCGeomServ::name(), RPCDetId::region(), RPCDetId::ring(), RPCDetId::sector(), DQMStore::setCurrentFolder(), RPCDetId::station(), and stripsInRoll().

Referenced by analyze().

00010                                                                                                       {
00011   map<string, MonitorElement*> meMap;  
00012 
00013   string ringType;
00014   int ring;
00015   if(detId.region() == 0) {
00016       ringType = "Wheel";  
00017     ring = detId.ring();
00018   }else if (detId.region() == -1){  
00019     ringType =  "Disk";
00020     ring = detId.region()*detId.station();
00021   }else {
00022     ringType =  "Disk";
00023     ring = detId.station();
00024   }
00025 
00026   RPCBookFolderStructure *  folderStr = new RPCBookFolderStructure();
00027   string folder = "RPC/RecHits/" +  folderStr->folderStructure(detId);
00028 
00029   dbe->setCurrentFolder(folder);
00030   
00031   //get number of strips in current roll
00032   int nstrips = this->stripsInRoll(detId, iSetup);
00033   if (nstrips == 0 ) nstrips = 1;
00034 
00036   RPCGeomServ RPCname(detId);
00037   string nameRoll = RPCname.name();
00038  
00039   stringstream os;
00040   os.str("");
00041   os<<"Occupancy_"<<nameRoll;
00042   meMap[os.str()] = dbe->book1D(os.str(), os.str(), nstrips, 0.5, nstrips+0.5);
00043 
00044   if (dqmexpert) {    
00045     os.str("");
00046     os<<"BXN_"<<nameRoll;
00047     meMap[os.str()] = dbe->book1D(os.str(), os.str(), 21, -10.5, 10.5);
00048 
00049     os.str("");
00050     os<<"ClusterSize_"<<nameRoll;
00051     meMap[os.str()] = dbe->book1D(os.str(), os.str(), 20, 0.5, 20.5);
00052   
00053     os.str("");
00054     os<<"NumberOfClusters_"<<nameRoll;
00055     meMap[os.str()] = dbe->book1D(os.str(), os.str(), 10, 0.5, 10.5);
00056 
00057     os.str("");
00058     os<<"NumberOfDigi_"<<nameRoll;
00059     meMap[os.str()] = dbe->book1D(os.str(), os.str(), 10, 0.5, 10.5);
00060         
00061     os.str("");
00062     os<<"CrossTalkLow_"<<nameRoll;
00063     meMap[os.str()] = dbe->book1D(os.str(), os.str(),nstrips, 0.5,nstrips+0.5 );
00064     
00065     os.str("");
00066     os<<"CrossTalkHigh_"<<nameRoll;
00067     meMap[os.str()] = dbe->book1D(os.str(), os.str(), nstrips, 0.5, nstrips+0.5);
00068     
00069     os.str("");
00070     os<<"BXWithData_"<<nameRoll;
00071     meMap[os.str()] = dbe->book1D(os.str(), os.str(), 10, 0.5, 10.5);
00072 
00074     os.str("");
00075     os<<"RecHitXPosition_"<<nameRoll;
00076     meMap[os.str()] = dbe->book1D(os.str(), os.str(), 80, -120, 120);
00077     
00078     os.str("");
00079     os<<"RecHitDX_"<<nameRoll;
00080     meMap[os.str()] = dbe->book1D(os.str(), os.str(), 30, -10, 10);
00081     
00082     os.str("");
00083     os<<"RecHitCounter_"<<nameRoll;
00084     meMap[os.str()] = dbe->book1D(os.str(), os.str(),20,0.5,20.5);
00085   }
00086   
00087   if (dqmsuperexpert) {    
00088     os.str("");
00089     os<<"BXN_vs_strip_"<<nameRoll;
00090     meMap[os.str()] = dbe->book2D(os.str(), os.str(),  nstrips , 0.5, nstrips+0.5 , 21, -10.5, 10.5);
00091         
00092     os.str("");
00093     os<<"ClusterSize_vs_LowerSrip_"<<nameRoll;
00094     meMap[os.str()] = dbe->book2D(os.str(), os.str(),  nstrips, 0.5,  nstrips+0.5,11, 0.5, 11.5);
00095     
00096     os.str("");
00097     os<<"ClusterSize_vs_HigherStrip_"<<nameRoll;
00098     meMap[os.str()] = dbe->book2D(os.str(), os.str(), nstrips, 0.5,  nstrips+0.5,11, 0.5, 11.5);
00099     
00100     os.str("");
00101     os<<"ClusterSize_vs_Strip_"<<nameRoll;
00102     meMap[os.str()] = dbe->book2D(os.str(), os.str(),nstrips, 0.5, nstrips+0.5,11, 0.5, 11.5);
00103     
00104     os.str("");
00105     os<<"ClusterSize_vs_CentralStrip_"<<nameRoll;
00106     meMap[os.str()] = dbe->book2D(os.str(), os.str(), nstrips, 0.5, nstrips+0.5,11, 0.5, 11.5);
00107  
00109     os.str("");
00110     os<<"MissingHits_"<<nameRoll;
00111     meMap[os.str()] = dbe->book2D(os.str(), os.str(),nstrips , 0, nstrips, 2, 0.,2.);
00112     
00113     os.str("");
00114     os<<"RecHitX_vs_dx_"<<nameRoll;
00115     meMap[os.str()] = dbe->book2D(os.str(), os.str(),30, -100, 100,30,10,10);
00116   }
00117 
00118   MonitorElement * myMe;
00119 
00120   os.str("");
00121   if(detId.region()==0)
00122     os<<"RPC/RecHits/Barrel/Wheel_"<<ring<<"/SummaryBySectors/";
00123   else if (detId.region()==1)
00124     os<<"RPC/RecHits/Endcap+/Disk_"<<ring<<"/SummaryBySectors/";
00125   else 
00126     os<<"RPC/RecHits/Endcap-/Disk_"<<ring<<"/SummaryBySectors/";
00127   string WheelSummary = os.str();
00128   dbe->setCurrentFolder(WheelSummary);
00129   
00130   os.str("");
00131   os<<"Occupancy_"<<ringType<<"_"<<ring<<"_Sector_"<<detId.sector();
00132   myMe = dbe->get(WheelSummary+"/"+os.str());
00133   
00134   //check if ME for this sector have already been booked
00135   if(myMe)  meMap[os.str()]=myMe;
00136   else {
00137     if(detId.region()==0){
00138       if (detId.sector()==9 || detId.sector()==11)
00139         meMap[os.str()] = dbe->book2D(os.str(), os.str(), 96, 0.5,96.5, 15, 0.5, 15.5);
00140       else  if (detId.sector()==4) 
00141         meMap[os.str()] = dbe->book2D(os.str(), os.str(),  96, 0.5, 96.5, 21, 0.5, 21.5);
00142       else
00143         meMap[os.str()] = dbe->book2D(os.str(), os.str(), 96, 0.5,  96.5, 17, 0.5, 17.5);
00144     }else{//Endcap
00145         meMap[os.str()] = dbe->book2D(os.str(), os.str(), 32, 0.5, 32.5, 54, 0.5, 54.5);
00146     }
00147   }
00148   
00149   os.str("");
00150   os<<"BxDistribution_"<<ringType<<"_"<<ring<<"_Sector_"<<detId.sector();
00151   myMe = dbe->get(folder+"/"+os.str());
00152   if(myMe)  meMap[os.str()]=myMe;
00153   else meMap[os.str()] = dbe->book1D(os.str(), os.str(), 11, -5.5, 5.5);
00154 
00155   os.str("");
00156   os<<"NumberOfDigi_"<<ringType<<"_"<<ring<<"_Sector_"<<detId.sector();
00157   myMe = dbe->get(folder+"/"+os.str());
00158   if(myMe)  meMap[os.str()]=myMe;
00159   else meMap[os.str()] = dbe->book1D(os.str(), os.str(), 50, 0.5, 50.5);
00160 
00161   os.str("");
00162   os<<"ClusterSize_"<<ringType<<"_"<<ring<<"_Sector_"<<detId.sector();
00163   myMe = dbe->get(folder+"/"+os.str());
00164   if(myMe)  meMap[os.str()]=myMe;
00165   else meMap[os.str()] = dbe->book1D(os.str(), os.str(), 20, 0.5, 20.5);
00166 
00167   os.str("");
00168   os<<"NumberOfClusters_"<<ringType<<"_"<<ring<<"_Sector_"<<detId.sector();
00169   myMe = dbe->get(folder+"/"+os.str());
00170   if(myMe)  meMap[os.str()]=myMe;
00171   else  meMap[os.str()] = dbe->book1D(os.str(), os.str(), 50, 0.5, 50.5);
00172 
00173   os.str("");
00174   os<<"BXWithData_"<<ringType<<"_"<<ring<<"_Sector_"<<detId.sector();
00175   myMe = dbe->get(folder+"/"+os.str());
00176   if(myMe)  meMap[os.str()]=myMe;
00177   else  meMap[os.str()] = dbe->book1D(os.str(), os.str(), 10, 0.5, 10.5);
00178 
00179   return meMap;
00180 }

map< string, MonitorElement * > RPCMonitorDigi::bookRegionRing ( int  region,
int  ring 
)

Booking of MonitoringElemnt at Wheel/Disk level.

Definition at line 182 of file RPCBookDetUnitME.cc.

References DQMStore::book1D(), DQMStore::book2D(), dbe, GlobalHistogramsFolder, and DQMStore::setCurrentFolder().

Referenced by analyze().

00182                                                                                 {  
00183   map<string, MonitorElement*> meMap;  
00184   string ringType = (region ==  0)?"Wheel":"Disk";
00185 
00186   dbe->setCurrentFolder(GlobalHistogramsFolder);
00187   stringstream os;
00188 
00189   os<<"OccupancyXY_"<<ringType<<"_"<<ring;
00190   meMap[os.str()] = dbe->book2D(os.str(), os.str(),63, -800, 800, 63, -800, 800);
00191 
00192   os.str("");
00193   os<<"ClusterSize_"<<ringType<<"_"<<ring;
00194   meMap[os.str()] = dbe->book1D(os.str(), os.str(),20, 0.5, 20.5);
00195 
00196   os.str("");
00197   os<<"1DOccupancy_"<<ringType<<"_";
00198   if (region!=0){
00199     ring = ring*region;
00200     os<<ring;
00201     meMap[os.str()] = dbe->book1D(os.str(), os.str(), 6, 0.5, 6.5);
00202   }else{
00203      os<<ring;
00204      meMap[os.str()] = dbe->book1D(os.str(), os.str(), 12, 0.5, 12.5);
00205   }
00206 
00207   os.str("");
00208   os<<"BxDistribution_"<<ringType<<"_"<<ring;
00209   meMap[os.str()] = dbe->book1D(os.str(), os.str(), 11, -5.5, 5.5);
00210 
00211   return meMap; 
00212 }

void RPCMonitorDigi::endJob ( void   )  [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 149 of file RPCMonitorDigi.cc.

References dbe, RootFileName, DQMStore::save(), and saveRootFile.

00150 {
00151   if(saveRootFile) dbe->save(RootFileName); 
00152   dbe = 0;
00153 }

int RPCMonitorDigi::stripsInRoll ( RPCDetId id,
const edm::EventSetup iSetup 
) [private]

Definition at line 215 of file RPCBookDetUnitME.cc.

References edm::EventSetup::get(), and RPCRoll::nstrips().

Referenced by bookDetUnitME().

00215                                                                           {
00216   edm::ESHandle<RPCGeometry> rpcgeo;
00217   iSetup.get<MuonGeometryRecord>().get(rpcgeo);
00218 
00219   const RPCRoll * rpcRoll = rpcgeo->roll(id);
00220 
00221   if (rpcRoll)
00222     return  rpcRoll->nstrips();
00223   else 
00224     return 1;
00225 }


Member Data Documentation

MonitorElement* RPCMonitorDigi::BarrelOccupancy [private]

Definition at line 77 of file RPCMonitorDigi.h.

Referenced by analyze(), beginJob(), and beginRun().

MonitorElement* RPCMonitorDigi::ClusterSize_for_Barrel [private]

Definition at line 71 of file RPCMonitorDigi.h.

Referenced by analyze(), beginJob(), and beginRun().

MonitorElement* RPCMonitorDigi::ClusterSize_for_BarrelandEndcaps [private]

Definition at line 75 of file RPCMonitorDigi.h.

Referenced by analyze(), beginJob(), and beginRun().

MonitorElement* RPCMonitorDigi::ClusterSize_for_EndcapNegative [private]

Definition at line 73 of file RPCMonitorDigi.h.

Referenced by analyze(), beginJob(), and beginRun().

MonitorElement* RPCMonitorDigi::ClusterSize_for_EndcapPositive [private]

Definition at line 72 of file RPCMonitorDigi.h.

Referenced by analyze(), beginJob(), and beginRun().

int RPCMonitorDigi::counter [private]

Definition at line 56 of file RPCMonitorDigi.h.

Referenced by analyze().

DQMStore* RPCMonitorDigi::dbe [private]

DQM store.

Definition at line 58 of file RPCMonitorDigi.h.

Referenced by beginJob(), bookDetUnitME(), bookRegionRing(), and endJob().

bool RPCMonitorDigi::dqmexpert [private]

Definition at line 96 of file RPCMonitorDigi.h.

Referenced by analyze(), bookDetUnitME(), and RPCMonitorDigi().

bool RPCMonitorDigi::dqmshifter [private]

Definition at line 95 of file RPCMonitorDigi.h.

Referenced by RPCMonitorDigi().

bool RPCMonitorDigi::dqmsuperexpert [private]

Definition at line 97 of file RPCMonitorDigi.h.

Referenced by analyze(), bookDetUnitME(), and RPCMonitorDigi().

MonitorElement* RPCMonitorDigi::EndcapNegativeOccupancy [private]

Definition at line 79 of file RPCMonitorDigi.h.

Referenced by analyze(), beginJob(), and beginRun().

MonitorElement* RPCMonitorDigi::EndcapPositiveOccupancy [private]

Definition at line 78 of file RPCMonitorDigi.h.

Referenced by analyze(), beginJob(), and beginRun().

std::map<uint32_t,bool> RPCMonitorDigi::foundHitsInChamber [private]

Definition at line 99 of file RPCMonitorDigi.h.

Referenced by analyze(), and RPCMonitorDigi().

std::string RPCMonitorDigi::GlobalHistogramsFolder [private]

Definition at line 98 of file RPCMonitorDigi.h.

Referenced by beginJob(), and bookRegionRing().

std::map<uint32_t, std::map<std::string, MonitorElement*> > RPCMonitorDigi::meCollection [private]

Definition at line 85 of file RPCMonitorDigi.h.

Referenced by analyze(), and beginRun().

bool RPCMonitorDigi::mergeRuns_ [private]

Definition at line 89 of file RPCMonitorDigi.h.

Referenced by beginRun(), and RPCMonitorDigi().

std::map<std::pair<int,int>, std::map<std::string, MonitorElement*> > RPCMonitorDigi::meWheelDisk [private]

Definition at line 86 of file RPCMonitorDigi.h.

Referenced by analyze(), and beginRun().

std::string RPCMonitorDigi::nameInLog [private]

Definition at line 91 of file RPCMonitorDigi.h.

Referenced by analyze(), beginJob(), beginRun(), and RPCMonitorDigi().

MonitorElement* RPCMonitorDigi::NumberOfClusters_for_Barrel [private]

Definition at line 67 of file RPCMonitorDigi.h.

Referenced by analyze(), beginJob(), and beginRun().

MonitorElement* RPCMonitorDigi::NumberOfClusters_for_EndcapNegative [private]

Definition at line 69 of file RPCMonitorDigi.h.

Referenced by analyze(), beginJob(), and beginRun().

MonitorElement* RPCMonitorDigi::NumberOfClusters_for_EndcapPositive [private]

Definition at line 68 of file RPCMonitorDigi.h.

Referenced by analyze(), beginJob(), and beginRun().

MonitorElement* RPCMonitorDigi::NumberOfDigis_for_Barrel [private]

Definition at line 63 of file RPCMonitorDigi.h.

Referenced by analyze(), beginJob(), and beginRun().

MonitorElement* RPCMonitorDigi::NumberOfDigis_for_EndcapNegative [private]

Definition at line 65 of file RPCMonitorDigi.h.

Referenced by analyze(), beginJob(), and beginRun().

MonitorElement* RPCMonitorDigi::NumberOfDigis_for_EndcapPositive [private]

Definition at line 64 of file RPCMonitorDigi.h.

Referenced by analyze(), beginJob(), and beginRun().

std::string RPCMonitorDigi::RootFileName [private]

Definition at line 94 of file RPCMonitorDigi.h.

Referenced by endJob(), and RPCMonitorDigi().

MonitorElement* RPCMonitorDigi::SameBxDigisMeBarrel_ [private]

Definition at line 81 of file RPCMonitorDigi.h.

Referenced by analyze(), beginJob(), and beginRun().

MonitorElement* RPCMonitorDigi::SameBxDigisMeEndcapNegative_ [private]

Definition at line 83 of file RPCMonitorDigi.h.

MonitorElement* RPCMonitorDigi::SameBxDigisMeEndcapPositive_ [private]

Definition at line 82 of file RPCMonitorDigi.h.

bool RPCMonitorDigi::saveRootFile [private]

Definition at line 92 of file RPCMonitorDigi.h.

Referenced by endJob(), and RPCMonitorDigi().

int RPCMonitorDigi::saveRootFileEventsInterval [private]

Definition at line 93 of file RPCMonitorDigi.h.

Referenced by RPCMonitorDigi().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:31:01 2009 for CMSSW by  doxygen 1.5.4