CMS 3D CMS Logo

DTLocalTriggerEfficiencyTest Class Reference

* DQM Test Client More...

#include <DQM/DTMonitorClient/src/DTLocalTriggerEfficiencyTest.h>

Inheritance diagram for DTLocalTriggerEfficiencyTest:

DTLocalTriggerBaseTest edm::EDAnalyzer

List of all members.

Public Member Functions

 DTLocalTriggerEfficiencyTest (const edm::ParameterSet &ps)
 Constructor.
virtual ~DTLocalTriggerEfficiencyTest ()
 Destructor.

Protected Member Functions

void beginJob (const edm::EventSetup &c)
 Begin Job.
void bookChambHistos (DTChamberId chambId, std::string htype)
 Book the new MEs (for each chamber).
void makeEfficiencyME (TH1D *numerator, TH1D *denominator, MonitorElement *result)
 Compute efficiency plots.
void makeEfficiencyME2D (TH2F *numerator, TH2F *denominator, MonitorElement *result)
 Compute 2D efficiency plots.
void runClientDiagnostic ()
 DQM Client Diagnostic.

Private Attributes

std::map< uint32_t, std::map
< std::string, MonitorElement * > > 
chambME


Detailed Description

* DQM Test Client

Date
2008/10/07 14:26:43
Revision
1.3
Author:
C. Battilana S. Marcellini - INFN Bologna

Definition at line 20 of file DTLocalTriggerEfficiencyTest.h.


Constructor & Destructor Documentation

DTLocalTriggerEfficiencyTest::DTLocalTriggerEfficiencyTest ( const edm::ParameterSet ps  ) 

Constructor.

Definition at line 36 of file DTLocalTriggerEfficiencyTest.cc.

References DTLocalTriggerBaseTest::baseFolder, and DTLocalTriggerBaseTest::setConfig().

00036                                                                                    {
00037 
00038   setConfig(ps,"DTLocalTriggerEfficiency");
00039   baseFolder = "DT/03-LocalTrigger/";
00040 
00041 }

DTLocalTriggerEfficiencyTest::~DTLocalTriggerEfficiencyTest (  )  [virtual]

Destructor.

Definition at line 44 of file DTLocalTriggerEfficiencyTest.cc.

00044                                                            {
00045   
00046 }


Member Function Documentation

void DTLocalTriggerEfficiencyTest::beginJob ( const edm::EventSetup c  )  [protected, virtual]

Begin Job.

Reimplemented from DTLocalTriggerBaseTest.

Definition at line 49 of file DTLocalTriggerEfficiencyTest.cc.

References DTLocalTriggerBaseTest::beginJob(), bookChambHistos(), DTLocalTriggerBaseTest::bookSectorHistos(), DTLocalTriggerBaseTest::bookWheelHistos(), edm::ParameterSet::getUntrackedParameter(), DTLocalTriggerBaseTest::hwSource, DTLocalTriggerBaseTest::hwSources, DTLocalTriggerBaseTest::parameters, DTLocalTriggerBaseTest::trigSource, and DTLocalTriggerBaseTest::trigSources.

00049                                                                  {
00050   
00051   DTLocalTriggerBaseTest::beginJob(c);
00052 
00053 
00054   vector<string>::const_iterator iTr   = trigSources.begin();
00055   vector<string>::const_iterator trEnd = trigSources.end();
00056   vector<string>::const_iterator iHw   = hwSources.begin();
00057   vector<string>::const_iterator hwEnd = hwSources.end();
00058 
00059 
00060   //Booking
00061   if(parameters.getUntrackedParameter<bool>("staticBooking", true)){
00062     for (; iTr != trEnd; ++iTr){
00063       trigSource = (*iTr);
00064       for (; iHw != hwEnd; ++iHw){
00065         hwSource = (*iHw);
00066         // Loop over the TriggerUnits
00067         for (int wh=-2; wh<=2; ++wh){
00068           for (int sect=1; sect<=12; ++sect){
00069             for (int stat=1; stat<=4; ++stat){
00070               DTChamberId chId(wh,stat,sect);
00071               bookChambHistos(chId,"TrigEffPosvsAnglePhi");
00072               bookChambHistos(chId,"TrigEffPosvsAngleHHHLPhi");
00073               bookChambHistos(chId,"TrigEffPosPhi");
00074               bookChambHistos(chId,"TrigEffPosHHHLPhi");
00075               bookChambHistos(chId,"TrigEffAnglePhi");
00076               bookChambHistos(chId,"TrigEffAngleHHHLPhi");
00077               bookChambHistos(chId,"TrigEffPosvsAngleTheta");
00078               bookChambHistos(chId,"TrigEffPosvsAngleHTheta");
00079               bookChambHistos(chId,"TrigEffPosTheta");
00080               bookChambHistos(chId,"TrigEffPosHTheta");
00081               bookChambHistos(chId,"TrigEffAngleTheta");
00082               bookChambHistos(chId,"TrigEffAngleHTheta");
00083             }
00084             bookSectorHistos(wh,sect,"","TrigEffPhi");  
00085             bookSectorHistos(wh,sect,"","TrigEffTheta");  
00086           }
00087           bookWheelHistos(wh,"","TrigEffPhi");  
00088           bookWheelHistos(wh,"","TrigEffHHHLPhi");  
00089           bookWheelHistos(wh,"","TrigEffTheta");  
00090           bookWheelHistos(wh,"","TrigEffHTheta");  
00091         }
00092       }
00093     }
00094   }
00095   
00096 }

void DTLocalTriggerEfficiencyTest::bookChambHistos ( DTChamberId  chambId,
std::string  htype 
) [protected]

Book the new MEs (for each chamber).

Referenced by beginJob(), and runClientDiagnostic().

void DTLocalTriggerEfficiencyTest::makeEfficiencyME ( TH1D *  numerator,
TH1D *  denominator,
MonitorElement result 
) [protected]

Compute efficiency plots.

Definition at line 237 of file DTLocalTriggerEfficiencyTest.cc.

References error, MonitorElement::getTH1F(), and funct::sqrt().

Referenced by runClientDiagnostic().

00237                                                                                                              {
00238   
00239   TH1F* efficiency = result->getTH1F();
00240   efficiency->Divide(numerator,denominator,1,1,"");
00241   
00242   int nbins = efficiency->GetNbinsX();
00243   for (int bin=1; bin<=nbins; ++bin){
00244     float error = 0;
00245     float bineff = efficiency->GetBinContent(bin);
00246 
00247     if (denominator->GetBinContent(bin)){
00248       error = sqrt(bineff*(1-bineff)/denominator->GetBinContent(bin));
00249     }
00250     else {
00251       error = 1;
00252       efficiency->SetBinContent(bin,1.);
00253     }
00254  
00255     efficiency->SetBinError(bin,error);
00256   }
00257 
00258 }

void DTLocalTriggerEfficiencyTest::makeEfficiencyME2D ( TH2F *  numerator,
TH2F *  denominator,
MonitorElement result 
) [protected]

Compute 2D efficiency plots.

Definition at line 261 of file DTLocalTriggerEfficiencyTest.cc.

References error, MonitorElement::getTH2F(), and funct::sqrt().

Referenced by runClientDiagnostic().

00261                                                                                                                {
00262   
00263   TH2F* efficiency = result->getTH2F();
00264   efficiency->Divide(numerator,denominator,1,1,"");
00265   
00266   int nbinsx = efficiency->GetNbinsX();
00267   int nbinsy = efficiency->GetNbinsY();
00268   for (int binx=1; binx<=nbinsx; ++binx){
00269     for (int biny=1; biny<=nbinsy; ++biny){
00270       float error = 0;
00271       float bineff = efficiency->GetBinContent(binx,biny);
00272 
00273       if (denominator->GetBinContent(binx,biny)){
00274         error = sqrt(bineff*(1-bineff)/denominator->GetBinContent(binx,biny));
00275       }
00276       else {
00277         error = 1;
00278         efficiency->SetBinContent(binx,biny,0.);
00279       }
00280  
00281       efficiency->SetBinError(binx,biny,error);
00282     }
00283   }
00284 
00285 }    

void DTLocalTriggerEfficiencyTest::runClientDiagnostic (  )  [protected, virtual]

DQM Client Diagnostic.

Implements DTLocalTriggerBaseTest.

Definition at line 99 of file DTLocalTriggerEfficiencyTest.cc.

References bookChambHistos(), DTLocalTriggerBaseTest::bookSectorHistos(), DTLocalTriggerBaseTest::bookWheelHistos(), chambME, DTLocalTriggerBaseTest::dbe, end, find(), DTLocalTriggerBaseTest::fullName(), DQMStore::get(), DTLocalTriggerBaseTest::getMEName(), DTLocalTriggerBaseTest::hwSource, DTLocalTriggerBaseTest::hwSources, makeEfficiencyME(), makeEfficiencyME2D(), DetId::rawId(), DTLocalTriggerBaseTest::secME, funct::sqrt(), DTLocalTriggerBaseTest::trigSource, DTLocalTriggerBaseTest::trigSources, and DTLocalTriggerBaseTest::whME.

00099                                                        {
00100 
00101   // Loop over Trig & Hw sources
00102   for (vector<string>::const_iterator iTr = trigSources.begin(); iTr != trigSources.end(); ++iTr){
00103     trigSource = (*iTr);
00104     for (vector<string>::const_iterator iHw = hwSources.begin(); iHw != hwSources.end(); ++iHw){
00105       hwSource = (*iHw);
00106       // Loop over the TriggerUnits
00107       for (int stat=1; stat<=4; ++stat){
00108         for (int wh=-2; wh<=2; ++wh){
00109           for (int sect=1; sect<=12; ++sect){
00110             DTChamberId chId(wh,stat,sect);
00111             int sector_id = (wh+3)+(sect-1)*5;
00112             uint32_t indexCh = chId.rawId();
00113 
00114             // Perform Efficiency analysis (Phi+Segments 2D)
00115             TH2F * TrackPosvsAngle            = getHisto<TH2F>(dbe->get(getMEName("TrackPosvsAngle","Segment", chId)));
00116             TH2F * TrackPosvsAngleandTrig     = getHisto<TH2F>(dbe->get(getMEName("TrackPosvsAngleandTrig","Segment", chId)));
00117             TH2F * TrackPosvsAngleandTrigHHHL = getHisto<TH2F>(dbe->get(getMEName("TrackPosvsAngleandTrigHHHL","Segment", chId)));
00118             
00119             if (TrackPosvsAngle && TrackPosvsAngleandTrig && TrackPosvsAngleandTrigHHHL && TrackPosvsAngle->GetEntries()>1) {
00120               
00121               if( chambME[indexCh].find(fullName("TrigEffAnglePhi")) == chambME[indexCh].end()){
00122                 bookChambHistos(chId,"TrigEffPosvsAnglePhi");
00123                 bookChambHistos(chId,"TrigEffPosvsAngleHHHLPhi");
00124                 bookChambHistos(chId,"TrigEffPosPhi");
00125                 bookChambHistos(chId,"TrigEffPosHHHLPhi");
00126                 bookChambHistos(chId,"TrigEffAnglePhi");
00127                 bookChambHistos(chId,"TrigEffAngleHHHLPhi");
00128               }
00129               if( secME[sector_id].find(fullName("TrigEffPhi")) == secME[sector_id].end() ){
00130                 bookSectorHistos(wh,sect,"","TrigEffPhi");  
00131               }
00132               if( whME[wh].find(fullName("TrigEffPhi")) == whME[wh].end() ){
00133                 bookWheelHistos(wh,"","TrigEffPhi");  
00134                 bookWheelHistos(wh,"","TrigEffHHHLPhi");  
00135               }
00136 
00137               std::map<std::string,MonitorElement*> *innerME = &(secME[sector_id]);
00138               TH1D* TrackPos               = TrackPosvsAngle->ProjectionY();
00139               TH1D* TrackAngle             = TrackPosvsAngle->ProjectionX();
00140               TH1D* TrackPosandTrig        = TrackPosvsAngleandTrig->ProjectionY();
00141               TH1D* TrackAngleandTrig      = TrackPosvsAngleandTrig->ProjectionX();
00142               TH1D* TrackPosandTrigHHHL    = TrackPosvsAngleandTrigHHHL->ProjectionY();
00143               TH1D* TrackAngleandTrigHHHL  = TrackPosvsAngleandTrigHHHL->ProjectionX();
00144               float binEff     = float(TrackPosandTrig->GetEntries())/TrackPos->GetEntries();
00145               float binEffHHHL = float(TrackPosandTrigHHHL->GetEntries())/TrackPos->GetEntries();
00146               float binErr     = sqrt(binEff*(1-binEff)/TrackPos->GetEntries());
00147               float binErrHHHL = sqrt(binEffHHHL*(1-binEffHHHL)/TrackPos->GetEntries());
00148           
00149               MonitorElement* globalEff = innerME->find(fullName("TrigEffPhi"))->second;
00150               globalEff->setBinContent(stat,binEff);
00151               globalEff->setBinError(stat,binErr);
00152 
00153               innerME = &(whME[wh]);
00154               globalEff = innerME->find(fullName("TrigEffPhi"))->second;
00155               globalEff->setBinContent(sect,stat,binEff);
00156               globalEff->setBinError(sect,stat,binErr);
00157               globalEff = innerME->find(fullName("TrigEffHHHLPhi"))->second;
00158               globalEff->setBinContent(sect,stat,binEffHHHL);
00159               globalEff->setBinError(sect,stat,binErrHHHL);
00160           
00161           
00162               innerME = &(chambME[indexCh]);
00163               makeEfficiencyME(TrackPosandTrig,TrackPos,innerME->find(fullName("TrigEffPosPhi"))->second);
00164               makeEfficiencyME(TrackPosandTrigHHHL,TrackPos,innerME->find(fullName("TrigEffPosHHHLPhi"))->second);
00165               makeEfficiencyME(TrackAngleandTrig,TrackAngle,innerME->find(fullName("TrigEffAnglePhi"))->second);
00166               makeEfficiencyME(TrackAngleandTrigHHHL,TrackAngle,innerME->find(fullName("TrigEffAngleHHHLPhi"))->second);
00167               makeEfficiencyME2D(TrackPosvsAngleandTrig,TrackPosvsAngle,innerME->find(fullName("TrigEffPosvsAnglePhi"))->second);
00168               makeEfficiencyME2D(TrackPosvsAngleandTrigHHHL,TrackPosvsAngle,innerME->find(fullName("TrigEffPosvsAngleHHHLPhi"))->second);
00169              
00170             }
00171         
00172             // Perform Efficiency analysis (Theta+Segments)
00173             TH2F * TrackThetaPosvsAngle            = getHisto<TH2F>(dbe->get(getMEName("TrackThetaPosvsAngle","Segment", chId)));
00174             TH2F * TrackThetaPosvsAngleandTrig     = getHisto<TH2F>(dbe->get(getMEName("TrackThetaPosvsAngleandTrig","Segment", chId)));
00175             TH2F * TrackThetaPosvsAngleandTrigH    = getHisto<TH2F>(dbe->get(getMEName("TrackThetaPosvsAngleandTrigH","Segment", chId)));
00176             
00177             if (TrackThetaPosvsAngle && TrackThetaPosvsAngleandTrig && TrackThetaPosvsAngleandTrigH && TrackThetaPosvsAngle->GetEntries()>1) {
00178               
00179               if( chambME[indexCh].find(fullName("TrigEffAngleTheta")) == chambME[indexCh].end()){
00180                 bookChambHistos(chId,"TrigEffPosvsAngleTheta");
00181                 bookChambHistos(chId,"TrigEffPosvsAngleHTheta");
00182                 bookChambHistos(chId,"TrigEffPosTheta");
00183                 bookChambHistos(chId,"TrigEffPosHTheta");
00184                 bookChambHistos(chId,"TrigEffAngleTheta");
00185                 bookChambHistos(chId,"TrigEffAngleHTheta");
00186               }
00187               if( secME[sector_id].find(fullName("TrigEffTheta")) == secME[sector_id].end() ){
00188                 bookSectorHistos(wh,sect,"","TrigEffTheta");  
00189               }
00190               if( whME[wh].find(fullName("TrigEffTheta")) == whME[wh].end() ){
00191                 bookWheelHistos(wh,"","TrigEffTheta");  
00192                 bookWheelHistos(wh,"","TrigEffHTheta");  
00193               }
00194 
00195               std::map<std::string,MonitorElement*> *innerME = &(secME[sector_id]);
00196               TH1D* TrackThetaPos               = TrackThetaPosvsAngle->ProjectionY();
00197               TH1D* TrackThetaAngle             = TrackThetaPosvsAngle->ProjectionX();
00198               TH1D* TrackThetaPosandTrig        = TrackThetaPosvsAngleandTrig->ProjectionY();
00199               TH1D* TrackThetaAngleandTrig      = TrackThetaPosvsAngleandTrig->ProjectionX();
00200               TH1D* TrackThetaPosandTrigH       = TrackThetaPosvsAngleandTrigH->ProjectionY();
00201               TH1D* TrackThetaAngleandTrigH     = TrackThetaPosvsAngleandTrigH->ProjectionX();
00202               float binEff  = float(TrackThetaPosandTrig->GetEntries())/TrackThetaPos->GetEntries();
00203               float binErr  = sqrt(binEff*(1-binEff)/TrackThetaPos->GetEntries());
00204               float binEffH = float(TrackThetaPosandTrigH->GetEntries())/TrackThetaPos->GetEntries();
00205               float binErrH = sqrt(binEffH*(1-binEffH)/TrackThetaPos->GetEntries());
00206           
00207               MonitorElement* globalEff = innerME->find(fullName("TrigEffTheta"))->second;
00208               globalEff->setBinContent(stat,binEff);
00209               globalEff->setBinError(stat,binErr);
00210 
00211               innerME = &(whME[wh]);
00212               globalEff = innerME->find(fullName("TrigEffTheta"))->second;
00213               globalEff->setBinContent(sect,stat,binEff);
00214               globalEff->setBinError(sect,stat,binErr);
00215               globalEff = innerME->find(fullName("TrigEffHTheta"))->second;
00216               globalEff->setBinContent(sect,stat,binEffH);
00217               globalEff->setBinError(sect,stat,binErrH);
00218           
00219               innerME = &(chambME[indexCh]);
00220               makeEfficiencyME(TrackThetaPosandTrig,TrackThetaPos,innerME->find(fullName("TrigEffPosTheta"))->second);
00221               makeEfficiencyME(TrackThetaPosandTrigH,TrackThetaPos,innerME->find(fullName("TrigEffPosHTheta"))->second);
00222               makeEfficiencyME(TrackThetaAngleandTrig,TrackThetaAngle,innerME->find(fullName("TrigEffAngleTheta"))->second);
00223               makeEfficiencyME(TrackThetaAngleandTrigH,TrackThetaAngle,innerME->find(fullName("TrigEffAngleHTheta"))->second);
00224               makeEfficiencyME2D(TrackThetaPosvsAngleandTrig,TrackThetaPosvsAngle,innerME->find(fullName("TrigEffPosvsAngleTheta"))->second);
00225               makeEfficiencyME2D(TrackThetaPosvsAngleandTrigH,TrackThetaPosvsAngle,innerME->find(fullName("TrigEffPosvsAngleHTheta"))->second);      
00226             }
00227 
00228           }
00229         }
00230       }
00231     }
00232   }     
00233 
00234 }


Member Data Documentation

std::map<uint32_t,std::map<std::string,MonitorElement*> > DTLocalTriggerEfficiencyTest::chambME [private]

Definition at line 51 of file DTLocalTriggerEfficiencyTest.h.

Referenced by runClientDiagnostic().


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