CMS 3D CMS Logo

Public Member Functions | Protected Member Functions | Private Attributes

DTTriggerEfficiencyTest Class Reference

#include <DTTriggerEfficiencyTest.h>

Inheritance diagram for DTTriggerEfficiencyTest:
DTLocalTriggerBaseTest edm::EDAnalyzer

List of all members.

Public Member Functions

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

Protected Member Functions

void beginJob ()
 BeginJob.
void beginRun (const edm::Run &r, const edm::EventSetup &c)
 BeginRun.
void bookChambHistos (DTChamberId chambId, std::string htype, std::string folder="")
 Book the new MEs (for each chamber)
void bookHistos (std::string hTag, std::string folder)
 Book the new MEs (global)
void bookWheelHistos (int wheel, std::string hTag, std::string folder)
 Book the new MEs (for each wheel)
std::string getMEName (std::string histoTag, std::string folder, int wh)
 Get the ME name (by wheel)
void makeEfficiencyME (TH2F *numerator, TH2F *denominator, MonitorElement *result2DWh)
 Compute 2D efficiency plots.
void makeEfficiencyME (TH2F *numerator, TH2F *denominator, MonitorElement *result2DWh, MonitorElement *result1DWh, MonitorElement *result1D)
 Compute 1D/2D efficiency plots.
void runClientDiagnostic ()
 DQM Client Diagnostic.

Private Attributes

std::map< uint32_t, std::map
< std::string, MonitorElement * > > 
chambME
bool detailedPlots
std::map< int, std::map
< std::string, MonitorElement * > > 
EffDistrPerWh
std::map< std::string,
MonitorElement * > 
globalEffDistr
DTTrigGeomUtilstrigGeomUtils

Detailed Description

* DQM Test Client

Date:
2011/10/28 08:12:44
Revision:
1.4
Author:
C. Battilana - CIEMAT

Definition at line 22 of file DTTriggerEfficiencyTest.h.


Constructor & Destructor Documentation

DTTriggerEfficiencyTest::DTTriggerEfficiencyTest ( const edm::ParameterSet ps)

Constructor.

Definition at line 39 of file DTTriggerEfficiencyTest.cc.

References edm::ParameterSet::getUntrackedParameter().

                                                                         {

  setConfig(ps,"DTTriggerEfficiency");
  baseFolderDCC = "DT/03-LocalTrigger-DCC/";
  baseFolderDDU = "DT/04-LocalTrigger-DDU/";
  detailedPlots = ps.getUntrackedParameter<bool>("detailedAnalysis",true);

}
DTTriggerEfficiencyTest::~DTTriggerEfficiencyTest ( ) [virtual]

Destructor.

Definition at line 49 of file DTTriggerEfficiencyTest.cc.

                                                 {

}

Member Function Documentation

void DTTriggerEfficiencyTest::beginJob ( void  ) [protected, virtual]

BeginJob.

Reimplemented from DTLocalTriggerBaseTest.

Definition at line 54 of file DTTriggerEfficiencyTest.cc.

References bk::beginJob().

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

BeginRun.

Reimplemented from DTLocalTriggerBaseTest.

Definition at line 61 of file DTTriggerEfficiencyTest.cc.

References DTLocalTriggerBaseTest::beginRun(), bookHistos(), and Parameters::parameters.

                                                                            {

  DTLocalTriggerBaseTest::beginRun(r,c);
  trigGeomUtils = new DTTrigGeomUtils(muonGeom);

  vector<string>::const_iterator iTr   = trigSources.begin();
  vector<string>::const_iterator trEnd = trigSources.end();
  vector<string>::const_iterator iHw   = hwSources.begin();
  vector<string>::const_iterator hwEnd = hwSources.end();


  //Booking
  if(parameters.getUntrackedParameter<bool>("staticBooking", true)){
    for (; iTr != trEnd; ++iTr){
      trigSource = (*iTr);
      for (; iHw != hwEnd; ++iHw){
        hwSource = (*iHw);
        // Loop over the TriggerUnits
        bookHistos("TrigEffPhi","");
        bookHistos("TrigEffCorrPhi","");
        for (int wh=-2; wh<=2; ++wh){
          if (detailedPlots) {
            for (int sect=1; sect<=12; ++sect){
              for (int stat=1; stat<=4; ++stat){
                DTChamberId chId(wh,stat,sect);
                bookChambHistos(chId,"TrigEffPosvsAnglePhi","Segment");
                bookChambHistos(chId,"TrigEffPosvsAngleCorrPhi","Segment");
              }
            }
          }
          bookWheelHistos(wh,"TrigEffPhi","");  
          bookWheelHistos(wh,"TrigEffCorrPhi","");  
        }
      }
    }
  }

}
void DTTriggerEfficiencyTest::bookChambHistos ( DTChamberId  chambId,
std::string  htype,
std::string  folder = "" 
) [protected]

Book the new MEs (for each chamber)

Definition at line 296 of file DTLocalTriggerEfficiencyTest.cc.

References DQMStore::book1D(), DQMStore::book2D(), DTLocalTriggerBaseTest::category(), DTLocalTriggerEfficiencyTest::chambME, DTLocalTriggerBaseTest::dbe, DTLocalTriggerBaseTest::fullName(), DTLocalTriggerBaseTest::hwSource, LogTrace, max(), min, pileupCalc::nbins, DTTrigGeomUtils::phiRange(), DetId::rawId(), DTChamberId::sector(), DQMStore::setCurrentFolder(), relativeConstraints::station, DTChamberId::station(), DTLocalTriggerBaseTest::testName, DTTrigGeomUtils::thetaRange(), DTLocalTriggerBaseTest::topFolder(), DTLocalTriggerEfficiencyTest::trigGeomUtils, and DTChamberId::wheel().

                                                                                    {
  
  stringstream wheel; wheel << chambId.wheel();
  stringstream station; station << chambId.station();   
  stringstream sector; sector << chambId.sector();

  string fullType  = fullName(htype);
  bool isDCC = hwSource=="DCC" ;
  string HistoName = fullType + "_W" + wheel.str() + "_Sec" + sector.str() + "_St" + station.str();

  dbe->setCurrentFolder(topFolder(isDCC) + "Wheel" + wheel.str() +
                        "/Sector" + sector.str() +
                        "/Station" + station.str() + "/Segment");

  LogTrace(category()) << "[" << testName << "Test]: booking " + topFolder(isDCC) + "Wheel" << wheel.str() 
                       <<"/Sector" << sector.str() << "/Station" << station.str() << "/Segment/" << HistoName;

  
  uint32_t indexChId = chambId.rawId();
  if (htype.find("TrigEffAnglePhi") == 0){
    chambME[indexChId][fullType] = dbe->book1D(HistoName.c_str(),"Trigger efficiency vs angle of incidence (Phi)",16,-40.,40.);
  }
  else if (htype.find("TrigEffAngleHHHLPhi") == 0){
    chambME[indexChId][fullType] = dbe->book1D(HistoName.c_str(),"Trigger efficiency (HH/HL) vs angle of incidence (Phi)",16,-40.,40.);
  }
  else if (htype.find("TrigEffAngleTheta") == 0){
    chambME[indexChId][fullType] = dbe->book1D(HistoName.c_str(),"Trigger efficiency vs angle of incidence (Theta)",16,-40.,40.);
  }
  else if (htype.find("TrigEffAngleHTheta") == 0){
    chambME[indexChId][fullType] = dbe->book1D(HistoName.c_str(),"Trigger efficiency (H) vs angle of incidence (Theta)",16,-40.,40.);
  }
  else if (htype.find("TrigEffPosPhi") == 0 ){
    float min,max;
    int nbins;
    trigGeomUtils->phiRange(chambId,min,max,nbins);
    chambME[indexChId][fullType] = dbe->book1D(HistoName.c_str(),"Trigger efficiency vs position (Phi)",nbins,min,max);
  }
  else if (htype.find("TrigEffPosvsAnglePhi") == 0 ){
    float min,max;
    int nbins;
    trigGeomUtils->phiRange(chambId,min,max,nbins);
    chambME[indexChId][fullType] = dbe->book2D(HistoName.c_str(),"Trigger efficiency position vs angle (Phi)",16,-40.,40.,nbins,min,max);
  }
  else if (htype.find("TrigEffPosvsAngleHHHLPhi") == 0 ){
    float min,max;
    int nbins;
    trigGeomUtils->phiRange(chambId,min,max,nbins);
    chambME[indexChId][fullType] = dbe->book2D(HistoName.c_str(),"Trigger efficiency (HH/HL) pos vs angle (Phi)",16,-40.,40.,nbins,min,max);
  }
  else if (htype.find("TrigEffPosHHHLPhi") == 0 ){
    float min,max;
    int nbins;
    trigGeomUtils->phiRange(chambId,min,max,nbins);
    chambME[indexChId][fullType] = dbe->book1D(HistoName.c_str(),"Trigger efficiency (HH/HL) vs position (Phi)",nbins,min,max);
  }
  else if (htype.find("TrigEffPosTheta") == 0){
    float min,max;
    int nbins;
    trigGeomUtils->thetaRange(chambId,min,max,nbins);
    chambME[indexChId][fullType] = dbe->book1D(HistoName.c_str(),"Trigger efficiency vs position (Theta)",nbins,min,max);
  }
  else if (htype.find("TrigEffPosHTheta") == 0){
    float min,max;
    int nbins;
    trigGeomUtils->thetaRange(chambId,min,max,nbins);
    chambME[indexChId][fullType] = dbe->book1D(HistoName.c_str(),"Trigger efficiency (H) vs position (Theta)",nbins,min,max);
  }
  else if (htype.find("TrigEffPosvsAngleTheta") == 0 ){
    float min,max;
    int nbins;
    trigGeomUtils->thetaRange(chambId,min,max,nbins);
    chambME[indexChId][fullType] = dbe->book2D(HistoName.c_str(),"Trigger efficiency pos vs angle (Theta)",16,-40.,40.,nbins,min,max);
  }
  else if (htype.find("TrigEffPosvsAngleHTheta") == 0 ){
    float min,max;
    int nbins;
    trigGeomUtils->thetaRange(chambId,min,max,nbins);
    chambME[indexChId][fullType] = dbe->book2D(HistoName.c_str(),"Trigger efficiency (H) pos vs angle (Theta)",16,-40.,40.,nbins,min,max);
  }

}
void DTTriggerEfficiencyTest::bookHistos ( std::string  hTag,
std::string  folder 
) [protected]

Book the new MEs (global)

void DTTriggerEfficiencyTest::bookWheelHistos ( int  wheel,
std::string  hTag,
std::string  folder 
) [protected]

Book the new MEs (for each wheel)

Reimplemented from DTLocalTriggerBaseTest.

std::string DTTriggerEfficiencyTest::getMEName ( std::string  histoTag,
std::string  folder,
int  wh 
) [protected]

Get the ME name (by wheel)

Reimplemented from DTLocalTriggerBaseTest.

void DTTriggerEfficiencyTest::makeEfficiencyME ( TH2F *  numerator,
TH2F *  denominator,
MonitorElement result2DWh 
) [protected]

Compute 2D efficiency plots.

Definition at line 202 of file DTTriggerEfficiencyTest.cc.

References postValidation_cfi::efficiency, error, MonitorElement::getTH2F(), and mathSSE::sqrt().

                                                                                                            {

  TH2F* efficiency = result2DWh->getTH2F();
  efficiency->Divide(numerator,denominator,1,1,"");

  int nbinsx = efficiency->GetNbinsX();
  int nbinsy = efficiency->GetNbinsY();
  for (int binx=1; binx<=nbinsx; ++binx){
    for (int biny=1; biny<=nbinsy; ++biny){
      float error = 0;
      float bineff = efficiency->GetBinContent(binx,biny);

      if (denominator->GetBinContent(binx,biny)){
        error = sqrt(bineff*(1-bineff)/denominator->GetBinContent(binx,biny));
      }
      else {
        error = 1;
        efficiency->SetBinContent(binx,biny,0.);
      }

      efficiency->SetBinError(binx,biny,error);
    }
  }

}    
void DTTriggerEfficiencyTest::makeEfficiencyME ( TH2F *  numerator,
TH2F *  denominator,
MonitorElement result2DWh,
MonitorElement result1DWh,
MonitorElement result1D 
) [protected]

Compute 1D/2D efficiency plots.

Definition at line 173 of file DTTriggerEfficiencyTest.cc.

References postValidation_cfi::efficiency, error, MonitorElement::Fill(), MonitorElement::getTH2F(), and mathSSE::sqrt().

                                                                                                                                                                  {

  TH2F* efficiency = result2DWh->getTH2F();
  efficiency->Divide(numerator,denominator,1,1,"");

  int nbinsx = efficiency->GetNbinsX();
  int nbinsy = efficiency->GetNbinsY();
  for (int binx=1; binx<=nbinsx; ++binx){
    for (int biny=1; biny<=nbinsy; ++biny){
      float error = 0;
      float bineff = efficiency->GetBinContent(binx,biny);

      result1DWh->Fill(bineff);
      result1D->Fill(bineff);

      if (denominator->GetBinContent(binx,biny)){
        error = sqrt(bineff*(1-bineff)/denominator->GetBinContent(binx,biny));
      }
      else {
        error = 1;
        efficiency->SetBinContent(binx,biny,0.);
      }

      efficiency->SetBinError(binx,biny,error);
    }
  }

}
void DTTriggerEfficiencyTest::runClientDiagnostic ( ) [protected, virtual]

DQM Client Diagnostic.

Implements DTLocalTriggerBaseTest.

Definition at line 101 of file DTTriggerEfficiencyTest.cc.

References bookHistos(), spr::find(), newFWLiteAna::fullName, DetId::rawId(), and edm::second().

                                                  {

  // Loop over Trig & Hw sources
  for (vector<string>::const_iterator iTr = trigSources.begin(); iTr != trigSources.end(); ++iTr){
    trigSource = (*iTr);
    for (vector<string>::const_iterator iHw = hwSources.begin(); iHw != hwSources.end(); ++iHw){
      hwSource = (*iHw);
      // Loop over the TriggerUnits
      if( globalEffDistr.find(fullName("TrigEffPhi")) == globalEffDistr.end() ){
        bookHistos("TrigEffPhi","");
        bookHistos("TrigEffCorrPhi","");
      }
      for (int wh=-2; wh<=2; ++wh){

        TH2F * TrigEffDenum   = getHisto<TH2F>(dbe->get(getMEName("TrigEffDenum","Task",wh)));
        TH2F * TrigEffNum     = getHisto<TH2F>(dbe->get(getMEName("TrigEffNum","Task",wh)));
        TH2F * TrigEffCorrNum = getHisto<TH2F>(dbe->get(getMEName("TrigEffCorrNum","Task",wh)));

        if (TrigEffDenum && TrigEffNum && TrigEffCorrNum && TrigEffDenum->GetEntries()>1) {

          if( whME[wh].find(fullName("TrigEffPhi")) == whME[wh].end() ){
            bookWheelHistos(wh,"TrigEffPhi","");  
            bookWheelHistos(wh,"TrigEffCorrPhi","");  
          }

          MonitorElement* Eff1DAll_TrigEffPhi = (&globalEffDistr)->find(fullName("TrigEffPhi"))->second;
          MonitorElement* Eff1DAll_TrigEffCorrPhi = (&globalEffDistr)->find(fullName("TrigEffCorrPhi"))->second;

          MonitorElement* Eff1DWh_TrigEffPhi = (&(EffDistrPerWh[wh]))->find(fullName("TrigEffPhi"))->second;
          MonitorElement* Eff1DWh_TrigEffCorrPhi = (&(EffDistrPerWh[wh]))->find(fullName("TrigEffCorrPhi"))->second;

          MonitorElement* Eff2DWh_TrigEffPhi = (&(whME[wh]))->find(fullName("TrigEffPhi"))->second;
          MonitorElement* Eff2DWh_TrigEffCorrPhi = (&(whME[wh]))->find(fullName("TrigEffCorrPhi"))->second;

          makeEfficiencyME(TrigEffNum,TrigEffDenum,Eff2DWh_TrigEffPhi,Eff1DWh_TrigEffPhi,Eff1DAll_TrigEffPhi);
          makeEfficiencyME(TrigEffCorrNum,TrigEffDenum,Eff2DWh_TrigEffCorrPhi,Eff1DWh_TrigEffCorrPhi,Eff1DAll_TrigEffCorrPhi);

        }

        if (detailedPlots) {
          for (int stat=1; stat<=4; ++stat){
            for (int sect=1; sect<=12; ++sect){
              DTChamberId chId(wh,stat,sect);
              uint32_t indexCh = chId.rawId();

              // Perform Efficiency analysis (Phi+Segments 2D)
              TH2F * TrackPosvsAngle        = getHisto<TH2F>(dbe->get(getMEName("TrackPosvsAngle","Segment", chId)));
              TH2F * TrackPosvsAngleAnyQual = getHisto<TH2F>(dbe->get(getMEName("TrackPosvsAngleAnyQual","Segment", chId)));
              TH2F * TrackPosvsAngleCorr    = getHisto<TH2F>(dbe->get(getMEName("TrackPosvsAngleCorr","Segment", chId)));

              if (TrackPosvsAngle && TrackPosvsAngleAnyQual && TrackPosvsAngleCorr && TrackPosvsAngle->GetEntries()>1) {

                if( chambME[indexCh].find(fullName("TrigEffAnglePhi")) == chambME[indexCh].end()){
                  bookChambHistos(chId,"TrigEffPosvsAnglePhi","Segment");
                  bookChambHistos(chId,"TrigEffPosvsAngleCorrPhi","Segment");
                }

                std::map<std::string,MonitorElement*> *innerME = &(chambME[indexCh]);
                makeEfficiencyME(TrackPosvsAngleAnyQual,TrackPosvsAngle,innerME->find(fullName("TrigEffPosvsAnglePhi"))->second);
                makeEfficiencyME(TrackPosvsAngleCorr,TrackPosvsAngle,innerME->find(fullName("TrigEffPosvsAngleCorrPhi"))->second);

              }
            }
          }
        }
      }

    }
  }     

}

Member Data Documentation

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

Definition at line 67 of file DTTriggerEfficiencyTest.h.

Definition at line 69 of file DTTriggerEfficiencyTest.h.

std::map<int,std::map<std::string,MonitorElement*> > DTTriggerEfficiencyTest::EffDistrPerWh [private]

Definition at line 66 of file DTTriggerEfficiencyTest.h.

std::map<std::string, MonitorElement*> DTTriggerEfficiencyTest::globalEffDistr [private]

Definition at line 65 of file DTTriggerEfficiencyTest.h.

Definition at line 68 of file DTTriggerEfficiencyTest.h.