CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

HeavyFlavorHarvesting Class Reference

Inheritance diagram for HeavyFlavorHarvesting:
edm::EDAnalyzer

List of all members.

Public Member Functions

virtual void analyze (const edm::Event &event, const edm::EventSetup &eventSetup)
virtual void endRun (const edm::Run &, const edm::EventSetup &)
 HeavyFlavorHarvesting (const edm::ParameterSet &pset)
virtual ~HeavyFlavorHarvesting ()

Private Member Functions

void calculateEfficiency (const ParameterSet &pset)
void calculateEfficiency1D (TH1 *num, TH1 *den, string name)
void calculateEfficiency2D (TH2F *num, TH2F *den, string name)

Private Attributes

DQMStoredqmStore
const VParameterSet efficiencies
string myDQMrootFolder

Detailed Description

Definition at line 26 of file HeavyFlavorHarvesting.cc.


Constructor & Destructor Documentation

HeavyFlavorHarvesting::HeavyFlavorHarvesting ( const edm::ParameterSet pset)

Definition at line 44 of file HeavyFlavorHarvesting.cc.

                                                                       :
  myDQMrootFolder( pset.getUntrackedParameter<string>("MyDQMrootFolder") ),
  efficiencies( pset.getUntrackedParameter<VParameterSet>("Efficiencies") )
{
}
HeavyFlavorHarvesting::~HeavyFlavorHarvesting ( ) [virtual]

Definition at line 188 of file HeavyFlavorHarvesting.cc.

                                             {
}

Member Function Documentation

virtual void HeavyFlavorHarvesting::analyze ( const edm::Event event,
const edm::EventSetup eventSetup 
) [inline, virtual]

Implements edm::EDAnalyzer.

Definition at line 33 of file HeavyFlavorHarvesting.cc.

{};
void HeavyFlavorHarvesting::calculateEfficiency ( const ParameterSet pset) [private]

Definition at line 61 of file HeavyFlavorHarvesting.cc.

References calculateEfficiency1D(), calculateEfficiency2D(), dqmStore, DQMStore::get(), MonitorElement::getTH1(), edm::ParameterSet::getUntrackedParameter(), LogDebug, myDQMrootFolder, and DQMStore::setCurrentFolder().

Referenced by endRun().

                                                                       {
//get hold of numerator and denominator histograms
  vector<string> numDenEffMEnames = pset.getUntrackedParameter<vector<string> >("NumDenEffMEnames");
  if(numDenEffMEnames.size()!=3){
    LogDebug("HLTriggerOfflineHeavyFlavor") << "NumDenEffMEnames must have three names"<<endl;
    return;
  }
  string denMEname = myDQMrootFolder+"/"+numDenEffMEnames[1];
  string numMEname = myDQMrootFolder+"/"+numDenEffMEnames[0];
  MonitorElement *denME = dqmStore->get(denMEname);
  MonitorElement *numME = dqmStore->get(numMEname);
  if(denME==0 || numME==0){
    LogDebug("HLTriggerOfflineHeavyFlavor") << "Could not find MEs: "<<denMEname<<" or "<<numMEname<<endl;
    return;
  }
  TH1 *den = denME->getTH1();
  TH1 *num = numME->getTH1();
  //check compatibility of the histograms  
  if( den->GetNbinsX() != num->GetNbinsX() || den->GetNbinsY() != num->GetNbinsY() ||  den->GetNbinsZ() != num->GetNbinsZ() ){
    LogDebug("HLTriggerOfflineHeavyFlavor") << "Monitoring elements "<<numMEname<<" and "<<denMEname<<"are incompatible"<<endl;
    return;
  }
  //figure out the directory and efficiency name  
  string effName = numDenEffMEnames[2];
  string effDir = myDQMrootFolder;
  string::size_type slashPos = effName.rfind('/');
  if ( string::npos != slashPos ) {
    effDir += "/"+effName.substr(0, slashPos);
    effName.erase(0, slashPos+1);
  }
  dqmStore->setCurrentFolder(effDir);
  //calculate the efficiencies
  int dimensions = num->GetDimension();
  if(dimensions==1){
    calculateEfficiency1D( num, den, effName );
  }else if(dimensions==2){
    calculateEfficiency2D( (TH2F*)num, (TH2F*)den, effName );
    TH1D* numX = ((TH2F*)num)->ProjectionX(); 
    TH1D* denX = ((TH2F*)den)->ProjectionX();
    calculateEfficiency1D( numX, denX, effName+"X" );
    delete numX;
    delete denX;
    TH1D* numY = ((TH2F*)num)->ProjectionY();
    TH1D* denY = ((TH2F*)den)->ProjectionY();
    calculateEfficiency1D( numY, denY, effName+"Y" );
    delete numY;
    delete denY;
  }else{
    return;
  }
}
void HeavyFlavorHarvesting::calculateEfficiency1D ( TH1 *  num,
TH1 *  den,
string  name 
) [private]

Definition at line 113 of file HeavyFlavorHarvesting.cc.

References DQMStore::bookProfile(), dqmStore, alignCSCRings::e, interpolateCardsSimple::eff, i, and mathSSE::sqrt().

Referenced by calculateEfficiency().

                                                                                     {
  TProfile* eff;
  if(num->GetXaxis()->GetXbins()->GetSize()==0){
    eff = new TProfile(effName.c_str(),effName.c_str(),num->GetXaxis()->GetNbins(),num->GetXaxis()->GetXmin(),num->GetXaxis()->GetXmax());
  }else{
    eff = new TProfile(effName.c_str(),effName.c_str(),num->GetXaxis()->GetNbins(),num->GetXaxis()->GetXbins()->GetArray());
  }
  eff->SetTitle(effName.c_str());
  eff->SetXTitle( num->GetXaxis()->GetTitle() );
  eff->SetYTitle("Efficiency");
  eff->SetOption("PE");
  eff->SetLineColor(2);
  eff->SetLineWidth(2);
  eff->SetMarkerStyle(20);
  eff->SetMarkerSize(0.8);
  eff->GetYaxis()->SetRangeUser(-0.001,1.001);
  eff->SetStats(kFALSE);
  for(int i=1;i<=num->GetNbinsX();i++){
    double e, low, high;
#if ROOT_VERSION_CODE >= ROOT_VERSION(5,27,0)
    if (int(den->GetBinContent(i))>0.) e= double(num->GetBinContent(i))/double(den->GetBinContent(i));
    else e=0.;
    low=TEfficiency::Wilson((double)den->GetBinContent(i),(double)num->GetBinContent(i),0.683,false);
    high=TEfficiency::Wilson((double)den->GetBinContent(i),(double)num->GetBinContent(i),0.683,true);
#else
    Efficiency( (double)num->GetBinContent(i), (double)den->GetBinContent(i), 0.683, e, low, high );
#endif
    double err = e-low>high-e ? e-low : high-e;
    //here is the trick to store info in TProfile:
    eff->SetBinContent( i, e );
    eff->SetBinEntries( i, 1 );
    eff->SetBinError( i, sqrt(e*e+err*err) );
  }
  dqmStore->bookProfile(effName,eff);
  delete eff;
}
void HeavyFlavorHarvesting::calculateEfficiency2D ( TH2F *  num,
TH2F *  den,
string  name 
) [private]

Definition at line 150 of file HeavyFlavorHarvesting.cc.

References DQMStore::bookProfile2D(), dqmStore, alignCSCRings::e, interpolateCardsSimple::eff, i, and mathSSE::sqrt().

Referenced by calculateEfficiency().

                                                                                       {
  TProfile2D* eff;
  if(num->GetXaxis()->GetXbins()->GetSize()==0 && num->GetYaxis()->GetXbins()->GetSize()==0){
    eff = new TProfile2D(effName.c_str(),effName.c_str(),num->GetXaxis()->GetNbins(),num->GetXaxis()->GetXmin(),num->GetXaxis()->GetXmax(),num->GetYaxis()->GetNbins(),num->GetYaxis()->GetXmin(),num->GetYaxis()->GetXmax());
  }else if(num->GetXaxis()->GetXbins()->GetSize()!=0 && num->GetYaxis()->GetXbins()->GetSize()==0){
    eff = new TProfile2D(effName.c_str(),effName.c_str(),num->GetXaxis()->GetNbins(),num->GetXaxis()->GetXbins()->GetArray(),num->GetYaxis()->GetNbins(),num->GetYaxis()->GetXmin(),num->GetYaxis()->GetXmax());
  }else if(num->GetXaxis()->GetXbins()->GetSize()==0 && num->GetYaxis()->GetXbins()->GetSize()!=0){
    eff = new TProfile2D(effName.c_str(),effName.c_str(),num->GetXaxis()->GetNbins(),num->GetXaxis()->GetXmin(),num->GetXaxis()->GetXmax(),num->GetYaxis()->GetNbins(),num->GetYaxis()->GetXbins()->GetArray());
  }else{
    eff = new TProfile2D(effName.c_str(),effName.c_str(),num->GetXaxis()->GetNbins(),num->GetXaxis()->GetXbins()->GetArray(),num->GetYaxis()->GetNbins(),num->GetYaxis()->GetXbins()->GetArray());
  }
  eff->SetTitle(effName.c_str());
  eff->SetXTitle( num->GetXaxis()->GetTitle() );
  eff->SetYTitle( num->GetYaxis()->GetTitle() );
  eff->SetZTitle("Efficiency");
  eff->SetOption("colztexte");
  eff->GetZaxis()->SetRangeUser(-0.001,1.001);
  eff->SetStats(kFALSE);
  for(int i=0;i<num->GetSize();i++){
    double e, low, high;
#if ROOT_VERSION_CODE >= ROOT_VERSION(5,27,0)
    if (int(den->GetBinContent(i))>0.) e= double(num->GetBinContent(i))/double(den->GetBinContent(i));
    else e=0.;
    low=TEfficiency::Wilson((double)den->GetBinContent(i),(double)num->GetBinContent(i),0.683,false);
    high=TEfficiency::Wilson((double)den->GetBinContent(i),(double)num->GetBinContent(i),0.683,true);
#else
    Efficiency( (double)num->GetBinContent(i), (double)den->GetBinContent(i), 0.683, e, low, high );
#endif
    double err = e-low>high-e ? e-low : high-e;
    //here is the trick to store info in TProfile:
    eff->SetBinContent( i, e );
    eff->SetBinEntries( i, 1 );
    eff->SetBinError( i, sqrt(e*e+err*err) );
  }
  dqmStore->bookProfile2D(effName,eff);
  delete eff;
}
void HeavyFlavorHarvesting::endRun ( const edm::Run ,
const edm::EventSetup  
) [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 50 of file HeavyFlavorHarvesting.cc.

References calculateEfficiency(), dqmStore, efficiencies, and cppFunctionSkipper::operator.

                                                                       {
  dqmStore = Service<DQMStore>().operator->();
  if( !dqmStore ){
    LogError("HLTriggerOfflineHeavyFlavor") << "Could not find DQMStore service\n";
    return;
  }
  for(VParameterSet::const_iterator pset = efficiencies.begin(); pset!=efficiencies.end(); pset++){
    calculateEfficiency(*pset);
  }
}

Member Data Documentation

Definition at line 41 of file HeavyFlavorHarvesting.cc.

Referenced by endRun().

Definition at line 40 of file HeavyFlavorHarvesting.cc.

Referenced by calculateEfficiency().