CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/DQMOffline/Alignment/src/MuonAlignmentSummary.cc

Go to the documentation of this file.
00001 /*
00002  *  DQM client for muon alignment summary
00003  *
00004  *  $Date: 2010/03/29 13:18:44 $
00005  *  $Revision: 1.6 $
00006  *  \author J. Fernandez - Univ. Oviedo <Javier.Fernandez@cern.ch>
00007  */
00008 
00009 #include "DQMOffline/Alignment/interface/MuonAlignmentSummary.h"
00010 
00011 
00012 
00013 
00014 MuonAlignmentSummary::MuonAlignmentSummary(const edm::ParameterSet& pSet) {
00015 
00016     parameters = pSet;
00017 
00018     meanPositionRange = parameters.getUntrackedParameter<double>("meanPositionRange");
00019     rmsPositionRange = parameters.getUntrackedParameter<double>("rmsPositionRange");
00020     meanAngleRange = parameters.getUntrackedParameter<double>("meanAngleRange");
00021     rmsAngleRange = parameters.getUntrackedParameter<double>("rmsAngleRange");
00022 
00023     doDT = parameters.getUntrackedParameter<bool>("doDT");
00024     doCSC = parameters.getUntrackedParameter<bool>("doCSC");
00025 
00026     MEFolderName = parameters.getParameter<std::string>("FolderName");
00027     topFolder << MEFolderName+"Alignment/Muon";
00028 
00029     if(!(doDT || doCSC) ) { 
00030         edm::LogError("MuonAlignmentSummary") <<" Error!! At least one Muon subsystem (DT or CSC) must be monitorized!!" << std::endl;
00031         edm::LogError("MuonAlignmentSummary") <<" Please enable doDT or doCSC to True in your python cfg file!!!" << std::endl;
00032         exit(1);
00033     }
00034   
00035     dbe = edm::Service<DQMStore>().operator->();
00036 
00037 }
00038 
00039 MuonAlignmentSummary::~MuonAlignmentSummary() { 
00040 }
00041 
00042 
00043 void MuonAlignmentSummary::beginRun(edm::Run const& run,edm::EventSetup const& iSetup) {
00044 
00045 
00046     metname = "MuonAlignmentSummary";
00047 
00048     LogTrace(metname)<<"[MuonAlignmentSummary] Parameters initialization";
00049   
00050 
00051         if (doDT){
00052         dbe->setCurrentFolder(topFolder.str()+"/DT");
00053         hLocalPositionDT=dbe->book2D("hLocalPositionDT","Local DT position (cm) absolute MEAN residuals;Sector;;cm", 14,1, 15,40,0,40);
00054         hLocalAngleDT=dbe->book2D("hLocalAngleDT","Local DT angle (rad) absolute MEAN residuals;Sector;;rad", 14,1, 15,40,0,40); 
00055         hLocalPositionRmsDT=dbe->book2D("hLocalPositionRmsDT","Local DT position (cm) RMS residuals;Sector;;cm", 14,1, 15,40,0,40);
00056         hLocalAngleRmsDT=dbe->book2D("hLocalAngleRmsDT","Local DT angle (rad) RMS residuals;Sector;;rad", 14,1, 15,40,0,40); 
00057 
00058         hLocalXMeanDT=dbe->book1D("hLocalXMeanDT","Distribution of absolute MEAN Local X (cm) residuals for DT;<X> (cm);number of chambers",100,0,meanPositionRange);
00059         hLocalXRmsDT=dbe->book1D("hLocalXRmsDT","Distribution of RMS Local X (cm) residuals for DT;X RMS (cm);number of chambers", 100,0,rmsPositionRange);
00060         hLocalYMeanDT=dbe->book1D("hLocalYMeanDT","Distribution of absolute MEAN Local Y (cm) residuals for DT;<Y> (cm);number of chambers", 100,0,meanPositionRange);
00061         hLocalYRmsDT=dbe->book1D("hLocalYRmsDT","Distribution of RMS Local Y (cm) residuals for DT;Y RMS (cm);number of chambers", 100,0,rmsPositionRange);
00062 
00063         hLocalPhiMeanDT=dbe->book1D("hLocalPhiMeanDT","Distribution of absolute MEAN #phi (rad) residuals for DT;<#phi>(rad);number of chambers", 100,0,meanAngleRange);
00064         hLocalPhiRmsDT=dbe->book1D("hLocalPhiRmsDT","Distribution of RMS #phi (rad) residuals for DT;#phi RMS (rad);number of chambers", 100,0,rmsAngleRange);
00065         hLocalThetaMeanDT=dbe->book1D("hLocalThetaMeanDT","Distribution of absolute MEAN #theta (rad) residuals for DT;<#theta>(rad);number of chambers", 100,0,meanAngleRange);
00066         hLocalThetaRmsDT=dbe->book1D("hLocalThetaRmsDT","Distribution of RMS #theta (rad) residuals for DT;#theta RMS (rad);number of chambers",100,0,rmsAngleRange);
00067 
00068         hLocalPositionDT->Reset();
00069         hLocalAngleDT->Reset();
00070         hLocalPositionRmsDT->Reset();
00071         hLocalAngleRmsDT->Reset();
00072         hLocalXMeanDT->Reset();
00073         hLocalXRmsDT->Reset();
00074         hLocalYMeanDT->Reset();
00075         hLocalYRmsDT->Reset();
00076         hLocalPhiMeanDT->Reset();
00077         hLocalPhiRmsDT->Reset();
00078         hLocalThetaMeanDT->Reset();
00079         hLocalThetaRmsDT->Reset();
00080 
00081         }
00082         
00083         if (doCSC){
00084         
00085         dbe->setCurrentFolder(topFolder.str()+"/CSC");
00086         hLocalPositionCSC=dbe->book2D("hLocalPositionCSC","Local CSC position (cm) absolute MEAN residuals;Sector;;cm",36,1,37,40,0,40);
00087         hLocalAngleCSC=dbe->book2D("hLocalAngleCSC","Local CSC angle (rad) absolute MEAN residuals;Sector;;rad", 36,1,37,40,0,40); 
00088         hLocalPositionRmsCSC=dbe->book2D("hLocalPositionRmsCSC","Local CSC position (cm) RMS residuals;Sector;;cm", 36,1,37,40,0,40);
00089         hLocalAngleRmsCSC=dbe->book2D("hLocalAngleRmsCSC","Local CSC angle (rad) RMS residuals;Sector;;rad", 36,1,37,40,0,40); 
00090         
00091         hLocalXMeanCSC=dbe->book1D("hLocalXMeanCSC","Distribution of absolute MEAN Local X (cm) residuals for CSC;<X> (cm);number of chambers",100,0,meanPositionRange);
00092         hLocalXRmsCSC=dbe->book1D("hLocalXRmsCSC","Distribution of RMS Local X (cm) residuals for CSC;X RMS (cm);number of chambers", 100,0,rmsPositionRange);
00093         hLocalYMeanCSC=dbe->book1D("hLocalYMeanCSC","Distribution of absolute MEAN Local Y (cm) residuals for CSC;<Y> (cm);number of chambers", 100,0,meanPositionRange);
00094         hLocalYRmsCSC=dbe->book1D("hLocalYRmsCSC","Distribution of RMS Local Y (cm) residuals for CSC;Y RMS (cm);number of chambers", 100,0,rmsPositionRange);
00095 
00096         hLocalPhiMeanCSC=dbe->book1D("hLocalPhiMeanCSC","Distribution of absolute MEAN #phi (rad) residuals for CSC;<#phi>(rad);number of chambers", 100,0,meanAngleRange);
00097         hLocalPhiRmsCSC=dbe->book1D("hLocalPhiRmsCSC","Distribution of RMS #phi (rad) residuals for CSC;#phi RMS (rad);number of chambers", 100,0,rmsAngleRange);
00098         hLocalThetaMeanCSC=dbe->book1D("hLocalThetaMeanCSC","Distribution of absolute MEAN #theta (rad) residuals for CSC;<#theta>(rad);number of chambers", 100,0,meanAngleRange);
00099         hLocalThetaRmsCSC=dbe->book1D("hLocalThetaRmsCSC","Distribution of RMS #theta (rad) residuals for CSC;#theta RMS (rad);number of chambers",100,0,rmsAngleRange);
00100 
00101         hLocalPositionCSC->Reset();
00102         hLocalAngleCSC->Reset();
00103         hLocalPositionRmsCSC->Reset();
00104         hLocalAngleRmsCSC->Reset();
00105         hLocalXMeanCSC->Reset();
00106         hLocalXRmsCSC->Reset();
00107         hLocalYMeanCSC->Reset();
00108         hLocalYRmsCSC->Reset();
00109         hLocalPhiMeanCSC->Reset();
00110         hLocalPhiRmsCSC->Reset();
00111         hLocalThetaMeanCSC->Reset();
00112         hLocalThetaRmsCSC->Reset();
00113         
00114 
00115         }
00116 }
00117 
00118 
00119 void MuonAlignmentSummary::endRun(edm::Run const& run, edm::EventSetup const& iSetup) {
00120 
00121 
00122     LogTrace(metname)<<"[MuonAlignmentSummary] Saving the histos";
00123 
00124     char binLabel[15];
00125 
00126     for (int station=-4; station<5;station++){
00127         if (doDT){
00128             if(station>0){
00129                                         
00130                 for(int wheel=-2;wheel<3;wheel++){
00131                         
00132                     for (int sector=1;sector<15;sector++){
00133                                 
00134                         if(!((sector==13 || sector ==14) && station!=4)){
00135                                         
00136                             std::stringstream Wheel; Wheel<<wheel;
00137                             std::stringstream Station; Station<<station;
00138                             std::stringstream Sector; Sector<<sector;
00139                                                                                 
00140                             std::string nameOfHistoLocalX="ResidualLocalX_W"+Wheel.str()+"MB"+Station.str()+"S"+Sector.str();
00141                             std::string nameOfHistoLocalPhi= "ResidualLocalPhi_W"+Wheel.str()+"MB"+Station.str()+"S"+Sector.str();
00142                             std::string nameOfHistoLocalTheta= "ResidualLocalTheta_W"+Wheel.str()+"MB"+Station.str()+"S"+Sector.str();
00143                             std::string nameOfHistoLocalY= "ResidualLocalY_W"+Wheel.str()+"MB"+Station.str()+"S"+Sector.str();
00144                                           
00145                             std::string path= topFolder.str()+
00146                                 "/DT/Wheel"+Wheel.str()+
00147                                 "/Station"+Station.str()+
00148                                 "/Sector"+Sector.str()+"/";
00149                                         
00150                             std::string histo = path+nameOfHistoLocalX;
00151 
00152                                         Int_t nstation=station - 1;
00153                                         Int_t nwheel=wheel+2;
00154                             MonitorElement * localX =dbe->get(histo);
00155                             if(localX){
00156                                         
00157                                 Double_t Mean=localX->getMean();
00158                                 Double_t Error=localX->getMeanError();
00159 
00160                                 Int_t ybin=1+nwheel*8+nstation*2;
00161                                 hLocalPositionDT->setBinContent(sector,ybin,fabs(Mean));
00162                                 sprintf(binLabel, "MB%d/%d_X",wheel, station );
00163                                 hLocalPositionDT->setBinLabel(ybin,binLabel,2);
00164                                 hLocalPositionRmsDT->setBinContent(sector,ybin,Error);
00165                                 hLocalPositionRmsDT->setBinLabel(ybin,binLabel,2);
00166 
00167                                 if(localX->getEntries()!=0){
00168                                     hLocalXMeanDT->Fill(fabs(Mean));
00169                                     hLocalXRmsDT->Fill(Error);}
00170                             }
00171 
00172                             histo = path+nameOfHistoLocalPhi;
00173                             MonitorElement * localPhi = dbe->get(histo);
00174                             if(localPhi){
00175  
00176                                 Double_t Mean=localPhi->getMean();
00177                                 Double_t Error=localPhi->getMeanError();
00178 
00179                                 Int_t ybin=1+nwheel*8+nstation*2;
00180                                 hLocalAngleDT->setBinContent(sector,ybin,fabs(Mean));
00181                                 sprintf(binLabel, "MB%d/%d_#phi", wheel,station );
00182                                 hLocalAngleDT->setBinLabel(ybin,binLabel,2);
00183                                 hLocalAngleRmsDT->setBinContent(sector,ybin,Error);
00184                                 hLocalAngleRmsDT->setBinLabel(ybin,binLabel,2);
00185                                         
00186                                 if(localPhi->getEntries()!=0){
00187                                     hLocalPhiMeanDT->Fill(fabs(Mean));
00188                                     hLocalPhiRmsDT->Fill(Error);}
00189                             }
00190 
00191                             if(station!=4){
00192 
00193                                 histo=path+nameOfHistoLocalY;
00194                                 MonitorElement * localY = dbe->get(histo);
00195                                 if(localY){
00196 
00197                                     Double_t Mean=localY->getMean();
00198                                     Double_t Error=localY->getMeanError();
00199 
00200                                     Int_t ybin=2+nwheel*8+nstation*2;
00201                                     hLocalPositionDT->setBinContent(sector,ybin,fabs(Mean));
00202                                     sprintf(binLabel, "MB%d/%d_Y", wheel,station );
00203                                     hLocalPositionDT->setBinLabel(ybin,binLabel,2);
00204                                     hLocalPositionRmsDT->setBinContent(sector,ybin,Error);
00205                                     hLocalPositionRmsDT->setBinLabel(ybin,binLabel,2);
00206                                     if(localY->getEntries()!=0){
00207                                         hLocalYMeanDT->Fill(fabs(Mean));
00208                                         hLocalYRmsDT->Fill(Error);}
00209                                 }
00210                                 histo = path+nameOfHistoLocalTheta;
00211                                 MonitorElement * localTheta = dbe->get(histo);
00212                                 if(localTheta){
00213                                     Double_t Mean=localTheta->getMean();
00214                                     Double_t Error=localTheta->getMeanError();
00215 
00216                                     Int_t ybin=2+nwheel*8+nstation*2;
00217                                     hLocalAngleDT->setBinContent(sector,ybin,fabs(Mean));
00218                                     sprintf(binLabel, "MB%d/%d_#theta",wheel,station );
00219                                     hLocalAngleDT->setBinLabel(ybin,binLabel,2);
00220                                     hLocalAngleRmsDT->setBinContent(sector,ybin,Error);
00221                                     hLocalAngleRmsDT->setBinLabel(ybin,binLabel,2);
00222                                     if(localTheta->getEntries()!=0){
00223                                         hLocalThetaMeanDT->Fill(fabs(Mean));
00224                                         hLocalThetaRmsDT->Fill(Error);}
00225                                 }
00226 
00227                             }// station != 4
00228 
00229                                         } //avoid non existing sectors
00230                     } //sector 
00231                 } //wheel
00232             } //station>0
00233         }// doDT
00234         
00235         if (doCSC){
00236             if(station!=0){
00237 
00238                 for(int ring=1;ring<5;ring++){
00239 
00240                     for(int chamber=1;chamber<37;chamber++){
00241                                 
00242                         if( !( ((abs(station)==2 || abs(station)==3 || abs(station)==4) && ring==1 && chamber>18) || 
00243                                 ((abs(station)==2 || abs(station)==3 || abs(station)==4) && ring>2)) ){
00244                             std::stringstream Ring; Ring<<ring;
00245                             std::stringstream Station; Station<<station;
00246                             std::stringstream Chamber; Chamber<<chamber;
00247                                                                                
00248                             std::string nameOfHistoLocalX="ResidualLocalX_ME"+Station.str()+"R"+Ring.str()+"C"+Chamber.str();
00249                             std::string nameOfHistoLocalPhi= "ResidualLocalPhi_ME"+Station.str()+"R"+Ring.str()+"C"+Chamber.str();
00250                             std::string nameOfHistoLocalTheta= "ResidualLocalTheta_ME"+Station.str()+"R"+Ring.str()+"C"+Chamber.str();
00251                             std::string nameOfHistoLocalY= "ResidualLocalY_ME"+Station.str()+"R"+Ring.str()+"C"+Chamber.str();
00252 
00253                             std::string path = topFolder.str()+
00254                                 "/CSC/Station"+Station.str()+
00255                                 "/Ring"+Ring.str()+
00256                                 "/Chamber"+Chamber.str()+"/";
00257 
00258                             Int_t ybin=abs(station)*2+ring; 
00259                             if(abs(station)==1) ybin=ring;
00260                             if (station>0) ybin=ybin+10;
00261                             else ybin = 11 -ybin;
00262                             std::string histo = path +  nameOfHistoLocalX;                      
00263                             MonitorElement * localX = dbe->get(histo);
00264                             if(localX){
00265 
00266                                 Double_t Mean=localX->getMean();
00267                                 Double_t Error=localX->getMeanError();
00268 
00269                                 Int_t ybin2=2*ybin-1;
00270                                 hLocalPositionCSC->setBinContent(chamber,ybin2,fabs(Mean));
00271                                 sprintf(binLabel, "ME%d/%d_X", station,ring );
00272                                 hLocalPositionCSC->setBinLabel(ybin2,binLabel,2);
00273                                 hLocalPositionRmsCSC->setBinContent(chamber,ybin2,Error);
00274                                 hLocalPositionRmsCSC->setBinLabel(ybin2,binLabel,2);
00275                                 if(localX->getEntries()!=0){
00276                                     hLocalXMeanCSC->Fill(fabs(Mean));
00277                                     hLocalXRmsCSC->Fill(Error);}
00278                             }
00279                             histo = path +      nameOfHistoLocalPhi;                    
00280 
00281                             MonitorElement * localPhi = dbe->get(histo);
00282                             if(localPhi){
00283 
00284 
00285                                 Double_t Mean=localPhi->getMean();
00286                                 Double_t Error=localPhi->getMeanError();
00287 
00288                                 Int_t ybin2=2*ybin-1;
00289                                 hLocalAngleCSC->setBinContent(chamber,ybin2,fabs(Mean));
00290                                 sprintf(binLabel, "ME%d/%d_#phi", station,ring );
00291                                 hLocalAngleCSC->setBinLabel(ybin2,binLabel,2);
00292                                 hLocalAngleRmsCSC->setBinContent(chamber,ybin2,Error);
00293                                 hLocalAngleRmsCSC->setBinLabel(ybin2,binLabel,2);
00294                                 if(localPhi->getEntries()!=0){
00295                                     hLocalPhiMeanCSC->Fill(fabs(Mean));
00296                                     hLocalPhiRmsCSC->Fill(Error);}
00297                             }
00298                             histo = path +      nameOfHistoLocalTheta;
00299                             MonitorElement * localTheta = dbe->get(histo);
00300                             if(localTheta){
00301 
00302 
00303                                 Double_t Mean=localTheta->getMean();
00304                                 Double_t Error=localTheta->getMeanError();
00305 
00306                                 Int_t ybin2=2*ybin;
00307                                 hLocalAngleCSC->setBinContent(chamber,ybin2,fabs(Mean));
00308                                 sprintf(binLabel, "ME%d/%d_#theta", station,ring );
00309                                 hLocalAngleCSC->setBinLabel(ybin2,binLabel,2);
00310                                 hLocalAngleRmsCSC->setBinContent(chamber,ybin2,Error);
00311                                 hLocalAngleRmsCSC->setBinLabel(ybin2,binLabel,2);
00312                                 if(localTheta->getEntries()!=0){
00313                                     hLocalThetaMeanCSC->Fill(fabs(Mean));
00314                                     hLocalThetaRmsCSC->Fill(Error);}
00315 
00316                             }                                   
00317                             histo = path +      nameOfHistoLocalY;
00318 
00319                             MonitorElement * localY = dbe->get(histo);
00320                             if(localY){
00321 
00322                                 Double_t Mean=localY->getMean();
00323                                 Double_t Error=localY->getMeanError();
00324 
00325                                 Int_t ybin2=2*ybin;
00326                                 hLocalPositionCSC->setBinContent(chamber,ybin2,fabs(Mean));
00327                                 sprintf(binLabel, "ME%d/%d_Y", station,ring );
00328                                 hLocalPositionCSC->setBinLabel(ybin2,binLabel,2);
00329                                 hLocalPositionRmsCSC->setBinContent(chamber,ybin2,Error);
00330                                 hLocalPositionRmsCSC->setBinLabel(ybin2,binLabel,2);
00331                                 if(localY->getEntries()!=0){
00332                                     hLocalYMeanCSC->Fill(fabs(Mean));
00333                                     hLocalYRmsCSC->Fill(Error);}
00334                             }
00335                         } //avoid non existing rings
00336                     } //chamber
00337                 } //ring
00338             } // station!=0
00339         }// doCSC
00340 
00341         } // loop on stations
00342 
00343 }
00344