#include <Histograms.h>
Public Member Functions | |
virtual void | Clear () |
virtual void | Fill (const reco::Particle::LorentzVector &p4, const double &resValue, const int charge) |
HFunctionResolution (TFile *outputFile, const TString &name, const double &ptMax=100, const int totBinsY=200) | |
virtual void | Write () |
~HFunctionResolution () | |
Protected Member Functions | |
int | getXindex (const double &x) const |
int | getYindex (const double &y) const |
Protected Attributes | |
double | deltaX_ |
double | deltaY_ |
TH1F * | hReso_ |
TProfile * | hResoVSEta_prof_ |
TProfile * | hResoVSPhi_prof_ |
TProfile * | hResoVSPhiMinus_prof_ |
TProfile * | hResoVSPhiPlus_prof_ |
TProfile * | hResoVSPt_Bar_prof_ |
TProfile * | hResoVSPt_Endc_17_prof_ |
TProfile * | hResoVSPt_Endc_20_prof_ |
TProfile * | hResoVSPt_Endc_24_prof_ |
TProfile * | hResoVSPt_Ovlap_prof_ |
TProfile * | hResoVSPt_prof_ |
TH2F * | hResoVSPtEta_ |
int ** | resoCount_ |
double ** | resoVsPtEta_ |
int | totBinsX_ |
int | totBinsY_ |
double | xMin_ |
double | yMin_ |
This histogram class fills a TProfile with the resolution evaluated from the resolution functions for single muon quantities. The resolution functions are used by MuScleFit to evaluate the mass resolution, which is the value seen by minuit and through it, corrections are evaluated.
In the end we will compare the histograms filled by this class (from the resolution function, reflecting the parameters changes done by minuit) with those filled comparing recoMuons with genMuons (the real resolutions).
Definition at line 1189 of file Histograms.h.
HFunctionResolution::HFunctionResolution | ( | TFile * | outputFile, |
const TString & | name, | ||
const double & | ptMax = 100 , |
||
const int | totBinsY = 200 |
||
) | [inline] |
Definition at line 1192 of file Histograms.h.
References deltaX_, deltaY_, hReso_, hResoVSEta_prof_, hResoVSPhi_prof_, hResoVSPhiMinus_prof_, hResoVSPhiPlus_prof_, hResoVSPt_Bar_prof_, hResoVSPt_Endc_17_prof_, hResoVSPt_Endc_20_prof_, hResoVSPt_Endc_24_prof_, hResoVSPt_Ovlap_prof_, hResoVSPt_prof_, hResoVSPtEta_, i, j, AlCaRecoCosmics_cfg::name, Histograms::name_, ExpressReco_HICollisions_FallBack::ptMax, resoCount_, resoVsPtEta_, totBinsX_, totBinsY_, xMin_, and yMin_.
: Histograms(outputFile, name) { name_ = name; totBinsX_ = 200; totBinsY_ = totBinsY; xMin_ = 0.; yMin_ = -3.0; double xMax = ptMax; double yMax = 3.0; deltaX_ = xMax - xMin_; deltaY_ = yMax - yMin_; hReso_ = new TH1F( name+"_Reso", "resolution", 1000, 0, 1 ); hResoVSPtEta_ = new TH2F( name+"_ResoVSPtEta", "resolution vs pt and #eta", totBinsX_, xMin_, xMax, totBinsY_, yMin_, yMax ); // Create and initialize the resolution arrays resoVsPtEta_ = new double*[totBinsX_]; resoCount_ = new int*[totBinsX_]; for( int i=0; i<totBinsX_; ++i ) { resoVsPtEta_[i] = new double[totBinsY_]; resoCount_[i] = new int[totBinsY_]; for( int j=0; j<totBinsY_; ++j ) { resoVsPtEta_[i][j] = 0; resoCount_[i][j] = 0; } } hResoVSPt_prof_ = new TProfile( name+"_ResoVSPt_prof", "resolution VS pt", totBinsX_, xMin_, xMax, yMin_, yMax); hResoVSPt_Bar_prof_ = new TProfile( name+"_ResoVSPt_Bar_prof", "resolution VS pt Barrel", totBinsX_, xMin_, xMax, yMin_, yMax); hResoVSPt_Endc_17_prof_ = new TProfile( name+"_ResoVSPt_Endc_1.7_prof", "resolution VS pt Endcap (1.4<eta<1.7)", totBinsX_, xMin_, xMax, yMin_, yMax); hResoVSPt_Endc_20_prof_ = new TProfile( name+"_ResoVSPt_Endc_2.0_prof", "resolution VS pt Endcap (1.7<eta<2.0)", totBinsX_, xMin_, xMax, yMin_, yMax); hResoVSPt_Endc_24_prof_ = new TProfile( name+"_ResoVSPt_Endc_2.4_prof", "resolution VS pt Endcap (2.0<eta<2.4)", totBinsX_, xMin_, xMax, yMin_, yMax); hResoVSPt_Ovlap_prof_ = new TProfile( name+"_ResoVSPt_Ovlap_prof", "resolution VS pt Overlap", totBinsX_, xMin_, xMax, yMin_, yMax); hResoVSEta_prof_ = new TProfile( name+"_ResoVSEta_prof", "resolution VS eta", totBinsY_, yMin_, yMax, 0, 1); //hResoVSTheta_prof_ = new TProfile( name+"_ResoVSTheta_prof", "resolution VS theta", 30, 0, TMath::Pi(), 0, 1); hResoVSPhiPlus_prof_ = new TProfile( name+"_ResoVSPhiPlus_prof", "resolution VS phi mu+", 14, -3.2, 3.2, 0, 1); hResoVSPhiMinus_prof_ = new TProfile( name+"_ResoVSPhiMinus_prof", "resolution VS phi mu-", 14, -3.2, 3.2, 0, 1); hResoVSPhi_prof_ = new TProfile( name+"_ResoVSPhi_prof", "resolution VS phi", 14, -3.2, 3.2, -1, 1); }
HFunctionResolution::~HFunctionResolution | ( | ) | [inline] |
Definition at line 1227 of file Histograms.h.
References hReso_, hResoVSEta_prof_, hResoVSPhi_prof_, hResoVSPhiMinus_prof_, hResoVSPhiPlus_prof_, hResoVSPt_Bar_prof_, hResoVSPt_Endc_17_prof_, hResoVSPt_Endc_20_prof_, hResoVSPt_Endc_24_prof_, hResoVSPt_Ovlap_prof_, hResoVSPt_prof_, hResoVSPtEta_, i, resoCount_, resoVsPtEta_, and totBinsX_.
{ delete hReso_; delete hResoVSPtEta_; // Free the resolution arrays for( int i=0; i<totBinsX_; ++i ) { delete[] resoVsPtEta_[i]; delete[] resoCount_[i]; } delete[] resoVsPtEta_; delete[] resoCount_; // Free the profile histograms delete hResoVSPt_prof_; delete hResoVSPt_Bar_prof_; delete hResoVSPt_Endc_17_prof_; delete hResoVSPt_Endc_20_prof_; delete hResoVSPt_Endc_24_prof_; delete hResoVSPt_Ovlap_prof_; delete hResoVSEta_prof_; //delete hResoVSTheta_prof_; delete hResoVSPhiPlus_prof_; delete hResoVSPhiMinus_prof_; delete hResoVSPhi_prof_; }
virtual void HFunctionResolution::Clear | ( | ) | [inline, virtual] |
Implements Histograms.
Definition at line 1330 of file Histograms.h.
References hReso_, hResoVSEta_prof_, hResoVSPhi_prof_, hResoVSPhiMinus_prof_, hResoVSPhiPlus_prof_, hResoVSPt_Bar_prof_, hResoVSPt_Endc_17_prof_, hResoVSPt_Endc_20_prof_, hResoVSPt_Endc_24_prof_, hResoVSPt_Ovlap_prof_, hResoVSPt_prof_, and hResoVSPtEta_.
{ hReso_->Clear(); hResoVSPtEta_->Clear(); hResoVSPt_prof_->Clear(); hResoVSPt_Bar_prof_->Clear(); hResoVSPt_Endc_17_prof_->Clear(); hResoVSPt_Endc_20_prof_->Clear(); hResoVSPt_Endc_24_prof_->Clear(); hResoVSPt_Ovlap_prof_->Clear(); hResoVSEta_prof_->Clear(); //hResoVSTheta_prof_->Clear(); hResoVSPhiPlus_prof_->Clear(); hResoVSPhiMinus_prof_->Clear(); hResoVSPhi_prof_->Clear(); }
virtual void HFunctionResolution::Fill | ( | const reco::Particle::LorentzVector & | p4, |
const double & | resValue, | ||
const int | charge | ||
) | [inline, virtual] |
Reimplemented from Histograms.
Reimplemented in HFunctionResolutionVarianceCheck.
Definition at line 1250 of file Histograms.h.
References getXindex(), getYindex(), hReso_, hResoVSEta_prof_, hResoVSPhi_prof_, hResoVSPhiMinus_prof_, hResoVSPhiPlus_prof_, hResoVSPt_Bar_prof_, hResoVSPt_Endc_17_prof_, hResoVSPt_Endc_20_prof_, hResoVSPt_Endc_24_prof_, hResoVSPt_Ovlap_prof_, hResoVSPt_prof_, resoCount_, resoVsPtEta_, totBinsX_, and totBinsY_.
{ if( resValue != resValue ) return; hReso_->Fill(resValue); // Fill the arrays with the resolution value and count int xIndex = getXindex(p4.Pt()); int yIndex = getYindex(p4.Eta()); if ( 0 <= xIndex && xIndex < totBinsX_ && 0 <= yIndex && yIndex < totBinsY_ ) { resoVsPtEta_[xIndex][yIndex] += resValue; // ATTENTION: we count only for positive muons because we are considering di-muon resonances // and we use this counter to compute the mean in the end. The resoVsPtEta value is filled with each muon, // but they must be considered independently (analogous to a factor 2) so in the end we would have // to divide by N/2, that is why we only increase the counter for half the muons. // if( charge > 0 ) // No more. Changing it here influences also other uses of this class. The macro FunctionTerms.cc // multiplies the terms by the 2 factor. resoCount_[xIndex][yIndex] += 1; // hResoVSPtEta->Fill(p4.Pt(), p4.Eta(), resValue); hResoVSPt_prof_->Fill(p4.Pt(),resValue); if(fabs(p4.Eta())<=0.9) hResoVSPt_Bar_prof_->Fill(p4.Pt(),resValue); else if(fabs(p4.Eta())>0.9 && fabs(p4.Eta())<=1.4 ) hResoVSPt_Ovlap_prof_->Fill(p4.Pt(),resValue); else if(fabs(p4.Eta())>1.4 && fabs(p4.Eta())<=1.7 ) hResoVSPt_Endc_17_prof_->Fill(p4.Pt(),resValue); else if(fabs(p4.Eta())>1.7 && fabs(p4.Eta())<=2.0 ) hResoVSPt_Endc_20_prof_->Fill(p4.Pt(),resValue); else hResoVSPt_Endc_24_prof_->Fill(p4.Pt(),resValue); hResoVSEta_prof_->Fill(p4.Eta(),resValue); //hResoVSTheta_prof_->Fill(p4.Theta(),resValue); if(charge>0) hResoVSPhiPlus_prof_->Fill(p4.Phi(),resValue); else if(charge<0) hResoVSPhiMinus_prof_->Fill(p4.Phi(),resValue); hResoVSPhi_prof_->Fill(p4.Phi(),resValue); } }
int HFunctionResolution::getXindex | ( | const double & | x | ) | const [inline, protected] |
int HFunctionResolution::getYindex | ( | const double & | y | ) | const [inline, protected] |
virtual void HFunctionResolution::Write | ( | ) | [inline, virtual] |
Implements Histograms.
Reimplemented in HFunctionResolutionVarianceCheck.
Definition at line 1291 of file Histograms.h.
References svgfig::canvas(), Histograms::histoDir_, hReso_, hResoVSEta_prof_, hResoVSPhi_prof_, hResoVSPhiMinus_prof_, hResoVSPhiPlus_prof_, hResoVSPt_Bar_prof_, hResoVSPt_Endc_17_prof_, hResoVSPt_Endc_20_prof_, hResoVSPt_Endc_24_prof_, hResoVSPt_Ovlap_prof_, hResoVSPt_prof_, hResoVSPtEta_, i, j, MultiGaussianStateTransform::N, Histograms::outputFile_, resoCount_, resoVsPtEta_, totBinsX_, and totBinsY_.
{ if(histoDir_ != 0) histoDir_->cd(); hReso_->Write(); for( int i=0; i<totBinsX_; ++i ) { for( int j=0; j<totBinsY_; ++j ) { int N = resoCount_[i][j]; // Fill with the mean value if( N != 0 ) hResoVSPtEta_->SetBinContent( i+1, j+1, resoVsPtEta_[i][j]/N ); else hResoVSPtEta_->SetBinContent( i+1, j+1, 0 ); } } hResoVSPtEta_->Write(); hResoVSPt_prof_->Write(); hResoVSPt_Bar_prof_->Write(); hResoVSPt_Endc_17_prof_->Write(); hResoVSPt_Endc_20_prof_->Write(); hResoVSPt_Endc_24_prof_->Write(); hResoVSPt_Ovlap_prof_->Write(); hResoVSEta_prof_->Write(); //hResoVSTheta_prof_->Write(); hResoVSPhiMinus_prof_->Write(); hResoVSPhiPlus_prof_->Write(); hResoVSPhi_prof_->Write(); TCanvas canvas(TString(hResoVSPtEta_->GetName())+"_canvas", TString(hResoVSPtEta_->GetTitle())+" canvas", 1000, 800); canvas.Divide(2); canvas.cd(1); hResoVSPtEta_->Draw("lego"); canvas.cd(2); hResoVSPtEta_->Draw("surf5"); canvas.Write(); hResoVSPtEta_->Write(); outputFile_->cd(); }
double HFunctionResolution::deltaX_ [protected] |
Definition at line 1370 of file Histograms.h.
Referenced by getXindex(), and HFunctionResolution().
double HFunctionResolution::deltaY_ [protected] |
Definition at line 1370 of file Histograms.h.
Referenced by getYindex(), and HFunctionResolution().
TH1F* HFunctionResolution::hReso_ [protected] |
Definition at line 1353 of file Histograms.h.
Referenced by Clear(), Fill(), HFunctionResolution(), Write(), and ~HFunctionResolution().
TProfile* HFunctionResolution::hResoVSEta_prof_ [protected] |
Definition at line 1363 of file Histograms.h.
Referenced by Clear(), Fill(), HFunctionResolution(), Write(), and ~HFunctionResolution().
TProfile* HFunctionResolution::hResoVSPhi_prof_ [protected] |
Definition at line 1367 of file Histograms.h.
Referenced by Clear(), Fill(), HFunctionResolution(), Write(), and ~HFunctionResolution().
TProfile* HFunctionResolution::hResoVSPhiMinus_prof_ [protected] |
Definition at line 1365 of file Histograms.h.
Referenced by Clear(), Fill(), HFunctionResolution(), Write(), and ~HFunctionResolution().
TProfile* HFunctionResolution::hResoVSPhiPlus_prof_ [protected] |
Definition at line 1366 of file Histograms.h.
Referenced by Clear(), Fill(), HFunctionResolution(), Write(), and ~HFunctionResolution().
TProfile* HFunctionResolution::hResoVSPt_Bar_prof_ [protected] |
Definition at line 1358 of file Histograms.h.
Referenced by Clear(), Fill(), HFunctionResolution(), Write(), and ~HFunctionResolution().
TProfile* HFunctionResolution::hResoVSPt_Endc_17_prof_ [protected] |
Definition at line 1359 of file Histograms.h.
Referenced by Clear(), Fill(), HFunctionResolution(), Write(), and ~HFunctionResolution().
TProfile* HFunctionResolution::hResoVSPt_Endc_20_prof_ [protected] |
Definition at line 1360 of file Histograms.h.
Referenced by Clear(), Fill(), HFunctionResolution(), Write(), and ~HFunctionResolution().
TProfile* HFunctionResolution::hResoVSPt_Endc_24_prof_ [protected] |
Definition at line 1361 of file Histograms.h.
Referenced by Clear(), Fill(), HFunctionResolution(), Write(), and ~HFunctionResolution().
TProfile* HFunctionResolution::hResoVSPt_Ovlap_prof_ [protected] |
Definition at line 1362 of file Histograms.h.
Referenced by Clear(), Fill(), HFunctionResolution(), Write(), and ~HFunctionResolution().
TProfile* HFunctionResolution::hResoVSPt_prof_ [protected] |
Definition at line 1357 of file Histograms.h.
Referenced by Clear(), Fill(), HFunctionResolution(), Write(), and ~HFunctionResolution().
TH2F* HFunctionResolution::hResoVSPtEta_ [protected] |
Definition at line 1354 of file Histograms.h.
Referenced by Clear(), HFunctionResolution(), Write(), and ~HFunctionResolution().
int** HFunctionResolution::resoCount_ [protected] |
Definition at line 1356 of file Histograms.h.
Referenced by Fill(), HFunctionResolution(), Write(), and ~HFunctionResolution().
double** HFunctionResolution::resoVsPtEta_ [protected] |
Definition at line 1355 of file Histograms.h.
Referenced by Fill(), HFunctionResolution(), Write(), and ~HFunctionResolution().
int HFunctionResolution::totBinsX_ [protected] |
Definition at line 1368 of file Histograms.h.
Referenced by Fill(), HFunctionResolutionVarianceCheck::Fill(), getXindex(), HFunctionResolution(), HFunctionResolutionVarianceCheck::HFunctionResolutionVarianceCheck(), Write(), HFunctionResolutionVarianceCheck::Write(), ~HFunctionResolution(), and HFunctionResolutionVarianceCheck::~HFunctionResolutionVarianceCheck().
int HFunctionResolution::totBinsY_ [protected] |
Definition at line 1368 of file Histograms.h.
Referenced by Fill(), HFunctionResolutionVarianceCheck::Fill(), getYindex(), HFunctionResolution(), HFunctionResolutionVarianceCheck::HFunctionResolutionVarianceCheck(), Write(), HFunctionResolutionVarianceCheck::Write(), and HFunctionResolutionVarianceCheck::~HFunctionResolutionVarianceCheck().
double HFunctionResolution::xMin_ [protected] |
Definition at line 1369 of file Histograms.h.
Referenced by getXindex(), and HFunctionResolution().
double HFunctionResolution::yMin_ [protected] |
Definition at line 1369 of file Histograms.h.
Referenced by getYindex(), and HFunctionResolution().