CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Protected Member Functions | Protected Attributes
HFunctionResolution Class Reference

#include <Histograms.h>

Inheritance diagram for HFunctionResolution:
Histograms HFunctionResolutionVarianceCheck

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 ()
 
- Public Member Functions inherited from Histograms
void declareHistograms ()
 
virtual void Fill (const reco::Particle::LorentzVector &p1, const reco::Particle::LorentzVector &p2)
 
virtual void Fill (const reco::Particle::LorentzVector &p1, const reco::Particle::LorentzVector &p2, const int charge, const double &weight=1.)
 
virtual void Fill (const CLHEP::HepLorentzVector &momentum1, const CLHEP::HepLorentzVector &momentum2)
 
virtual void Fill (const CLHEP::HepLorentzVector &momentum1, const CLHEP::HepLorentzVector &momentum2, const int charge, const double &weight=1.)
 
virtual void Fill (const CLHEP::HepLorentzVector &p1, const reco::Particle::LorentzVector &p2)
 
virtual void Fill (const reco::Particle::LorentzVector &p4, const double &weight=1.)
 
virtual void Fill (const reco::Particle::LorentzVector &p4, const double &genValue, const double recValue, const int charge)
 
virtual void Fill (const CLHEP::HepLorentzVector &p, const double &likeValue)
 
virtual void Fill (const unsigned int number)
 
virtual void Fill (const reco::Particle::LorentzVector &recoP1, const int charge1, const reco::Particle::LorentzVector &genP1, const reco::Particle::LorentzVector &recoP2, const int charge2, const reco::Particle::LorentzVector &genP2, const double &recoMass, const double &genMass)
 
virtual void Fill (const reco::Particle::LorentzVector &recoP1, const int charge1, const reco::Particle::LorentzVector &recoP2, const int charge2, const double &recoMass, const double &genMass)
 
virtual void Fill (const reco::Particle::LorentzVector &recoP1, const reco::Particle::LorentzVector &genP1, const reco::Particle::LorentzVector &recoP2, const reco::Particle::LorentzVector &genP2)
 
virtual void Fill (const double &x, const double &y)
 
virtual void Fill (const double &x, const double &y, const double &a, const double &b)
 
void fillEventInfo (int proc, int strk, int ntrkr)
 
void fillRecHistograms (const RecTrack_t &r)
 
void fillSimHistograms (const SimTrack_t &s)
 
void fillVzeroHistograms (const RecVzero_t &r, int part)
 
virtual double Get (const reco::Particle::LorentzVector &recoP1, const TString &covarianceName)
 
virtual TString GetName ()
 
 Histograms ()
 
 Histograms (const TString &name)
 
 Histograms (TFile *outputFile, const TString &name)
 
 Histograms (const edm::ParameterSet &pset)
 
virtual void SetWeight (double weight)
 
void writeHistograms ()
 
virtual ~Histograms ()
 
 ~Histograms ()
 

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_
 
- Protected Attributes inherited from Histograms
TDirectory * histoDir_
 
TString name_
 
TFile * outputFile_
 
double theWeight_
 

Detailed Description

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.

Constructor & Destructor Documentation

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, mergeVDriftHistosByStation::name, Histograms::name_, jptDQMConfig_cff::ptMax, resoCount_, resoVsPtEta_, totBinsX_, totBinsY_, xMin_, and yMin_.

1192  : Histograms(outputFile, name) {
1193  name_ = name;
1194  totBinsX_ = 200;
1195  totBinsY_ = totBinsY;
1196  xMin_ = 0.;
1197  yMin_ = -3.0;
1198  double xMax = ptMax;
1199  double yMax = 3.0;
1200  deltaX_ = xMax - xMin_;
1201  deltaY_ = yMax - yMin_;
1202  hReso_ = new TH1F( name+"_Reso", "resolution", 1000, 0, 1 );
1203  hResoVSPtEta_ = new TH2F( name+"_ResoVSPtEta", "resolution vs pt and #eta", totBinsX_, xMin_, xMax, totBinsY_, yMin_, yMax );
1204  // Create and initialize the resolution arrays
1205  resoVsPtEta_ = new double*[totBinsX_];
1206  resoCount_ = new int*[totBinsX_];
1207  for( int i=0; i<totBinsX_; ++i ) {
1208  resoVsPtEta_[i] = new double[totBinsY_];
1209  resoCount_[i] = new int[totBinsY_];
1210  for( int j=0; j<totBinsY_; ++j ) {
1211  resoVsPtEta_[i][j] = 0;
1212  resoCount_[i][j] = 0;
1213  }
1214  }
1215  hResoVSPt_prof_ = new TProfile( name+"_ResoVSPt_prof", "resolution VS pt", totBinsX_, xMin_, xMax, yMin_, yMax);
1216  hResoVSPt_Bar_prof_ = new TProfile( name+"_ResoVSPt_Bar_prof", "resolution VS pt Barrel", totBinsX_, xMin_, xMax, yMin_, yMax);
1217  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);
1218  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);
1219  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);
1220  hResoVSPt_Ovlap_prof_ = new TProfile( name+"_ResoVSPt_Ovlap_prof", "resolution VS pt Overlap", totBinsX_, xMin_, xMax, yMin_, yMax);
1221  hResoVSEta_prof_ = new TProfile( name+"_ResoVSEta_prof", "resolution VS eta", totBinsY_, yMin_, yMax, 0, 1);
1222  //hResoVSTheta_prof_ = new TProfile( name+"_ResoVSTheta_prof", "resolution VS theta", 30, 0, TMath::Pi(), 0, 1);
1223  hResoVSPhiPlus_prof_ = new TProfile( name+"_ResoVSPhiPlus_prof", "resolution VS phi mu+", 14, -3.2, 3.2, 0, 1);
1224  hResoVSPhiMinus_prof_ = new TProfile( name+"_ResoVSPhiMinus_prof", "resolution VS phi mu-", 14, -3.2, 3.2, 0, 1);
1225  hResoVSPhi_prof_ = new TProfile( name+"_ResoVSPhi_prof", "resolution VS phi", 14, -3.2, 3.2, -1, 1);
1226  }
int i
Definition: DBlmapReader.cc:9
TString name_
Definition: Histograms.h:106
TProfile * hResoVSPt_Endc_20_prof_
Definition: Histograms.h:1360
TProfile * hResoVSPt_Endc_17_prof_
Definition: Histograms.h:1359
TProfile * hResoVSPt_Bar_prof_
Definition: Histograms.h:1358
TProfile * hResoVSPt_prof_
Definition: Histograms.h:1357
int j
Definition: DBlmapReader.cc:9
TProfile * hResoVSPt_Ovlap_prof_
Definition: Histograms.h:1362
TProfile * hResoVSPt_Endc_24_prof_
Definition: Histograms.h:1361
TProfile * hResoVSPhi_prof_
Definition: Histograms.h:1367
TProfile * hResoVSPhiMinus_prof_
Definition: Histograms.h:1365
TProfile * hResoVSEta_prof_
Definition: Histograms.h:1363
TProfile * hResoVSPhiPlus_prof_
Definition: Histograms.h:1366
double ** resoVsPtEta_
Definition: Histograms.h:1355
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_.

1227  {
1228  delete hReso_;
1229  delete hResoVSPtEta_;
1230  // Free the resolution arrays
1231  for( int i=0; i<totBinsX_; ++i ) {
1232  delete[] resoVsPtEta_[i];
1233  delete[] resoCount_[i];
1234  }
1235  delete[] resoVsPtEta_;
1236  delete[] resoCount_;
1237  // Free the profile histograms
1238  delete hResoVSPt_prof_;
1239  delete hResoVSPt_Bar_prof_;
1240  delete hResoVSPt_Endc_17_prof_;
1241  delete hResoVSPt_Endc_20_prof_;
1242  delete hResoVSPt_Endc_24_prof_;
1243  delete hResoVSPt_Ovlap_prof_;
1244  delete hResoVSEta_prof_;
1245  //delete hResoVSTheta_prof_;
1246  delete hResoVSPhiPlus_prof_;
1247  delete hResoVSPhiMinus_prof_;
1248  delete hResoVSPhi_prof_;
1249  }
int i
Definition: DBlmapReader.cc:9
TProfile * hResoVSPt_Endc_20_prof_
Definition: Histograms.h:1360
TProfile * hResoVSPt_Endc_17_prof_
Definition: Histograms.h:1359
TProfile * hResoVSPt_Bar_prof_
Definition: Histograms.h:1358
TProfile * hResoVSPt_prof_
Definition: Histograms.h:1357
TProfile * hResoVSPt_Ovlap_prof_
Definition: Histograms.h:1362
TProfile * hResoVSPt_Endc_24_prof_
Definition: Histograms.h:1361
TProfile * hResoVSPhi_prof_
Definition: Histograms.h:1367
TProfile * hResoVSPhiMinus_prof_
Definition: Histograms.h:1365
TProfile * hResoVSEta_prof_
Definition: Histograms.h:1363
TProfile * hResoVSPhiPlus_prof_
Definition: Histograms.h:1366
double ** resoVsPtEta_
Definition: Histograms.h:1355

Member Function Documentation

virtual void HFunctionResolution::Clear ( )
inlinevirtual

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_.

1330  {
1331  hReso_->Clear();
1332  hResoVSPtEta_->Clear();
1333  hResoVSPt_prof_->Clear();
1334  hResoVSPt_Bar_prof_->Clear();
1335  hResoVSPt_Endc_17_prof_->Clear();
1336  hResoVSPt_Endc_20_prof_->Clear();
1337  hResoVSPt_Endc_24_prof_->Clear();
1338  hResoVSPt_Ovlap_prof_->Clear();
1339  hResoVSEta_prof_->Clear();
1340  //hResoVSTheta_prof_->Clear();
1341  hResoVSPhiPlus_prof_->Clear();
1342  hResoVSPhiMinus_prof_->Clear();
1343  hResoVSPhi_prof_->Clear();
1344  }
TProfile * hResoVSPt_Endc_20_prof_
Definition: Histograms.h:1360
TProfile * hResoVSPt_Endc_17_prof_
Definition: Histograms.h:1359
TProfile * hResoVSPt_Bar_prof_
Definition: Histograms.h:1358
TProfile * hResoVSPt_prof_
Definition: Histograms.h:1357
TProfile * hResoVSPt_Ovlap_prof_
Definition: Histograms.h:1362
TProfile * hResoVSPt_Endc_24_prof_
Definition: Histograms.h:1361
TProfile * hResoVSPhi_prof_
Definition: Histograms.h:1367
TProfile * hResoVSPhiMinus_prof_
Definition: Histograms.h:1365
TProfile * hResoVSEta_prof_
Definition: Histograms.h:1363
TProfile * hResoVSPhiPlus_prof_
Definition: Histograms.h:1366
virtual void HFunctionResolution::Fill ( const reco::Particle::LorentzVector p4,
const double &  resValue,
const int  charge 
)
inlinevirtual

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_.

Referenced by HFunctionResolutionVarianceCheck::Fill().

1250  {
1251  if( resValue != resValue ) return;
1252  hReso_->Fill(resValue);
1253 
1254  // Fill the arrays with the resolution value and count
1255  int xIndex = getXindex(p4.Pt());
1256  int yIndex = getYindex(p4.Eta());
1257  if ( 0 <= xIndex && xIndex < totBinsX_ && 0 <= yIndex && yIndex < totBinsY_ ) {
1258  resoVsPtEta_[xIndex][yIndex] += resValue;
1259  // ATTENTION: we count only for positive muons because we are considering di-muon resonances
1260  // and we use this counter to compute the mean in the end. The resoVsPtEta value is filled with each muon,
1261  // but they must be considered independently (analogous to a factor 2) so in the end we would have
1262  // to divide by N/2, that is why we only increase the counter for half the muons.
1263  // if( charge > 0 )
1264  // No more. Changing it here influences also other uses of this class. The macro FunctionTerms.cc
1265  // multiplies the terms by the 2 factor.
1266 
1267  resoCount_[xIndex][yIndex] += 1;
1268 
1269  // hResoVSPtEta->Fill(p4.Pt(), p4.Eta(), resValue);
1270  hResoVSPt_prof_->Fill(p4.Pt(),resValue);
1271  if(fabs(p4.Eta())<=0.9)
1272  hResoVSPt_Bar_prof_->Fill(p4.Pt(),resValue);
1273  else if(fabs(p4.Eta())>0.9 && fabs(p4.Eta())<=1.4 )
1274  hResoVSPt_Ovlap_prof_->Fill(p4.Pt(),resValue);
1275  else if(fabs(p4.Eta())>1.4 && fabs(p4.Eta())<=1.7 )
1276  hResoVSPt_Endc_17_prof_->Fill(p4.Pt(),resValue);
1277  else if(fabs(p4.Eta())>1.7 && fabs(p4.Eta())<=2.0 )
1278  hResoVSPt_Endc_20_prof_->Fill(p4.Pt(),resValue);
1279  else
1280  hResoVSPt_Endc_24_prof_->Fill(p4.Pt(),resValue);
1281  hResoVSEta_prof_->Fill(p4.Eta(),resValue);
1282  //hResoVSTheta_prof_->Fill(p4.Theta(),resValue);
1283  if(charge>0)
1284  hResoVSPhiPlus_prof_->Fill(p4.Phi(),resValue);
1285  else if(charge<0)
1286  hResoVSPhiMinus_prof_->Fill(p4.Phi(),resValue);
1287  hResoVSPhi_prof_->Fill(p4.Phi(),resValue);
1288  }
1289  }
double charge(const std::vector< uint8_t > &Ampls)
int getYindex(const double &y) const
Definition: Histograms.h:1350
TProfile * hResoVSPt_Endc_20_prof_
Definition: Histograms.h:1360
TProfile * hResoVSPt_Endc_17_prof_
Definition: Histograms.h:1359
TProfile * hResoVSPt_Bar_prof_
Definition: Histograms.h:1358
double p4[4]
Definition: TauolaWrapper.h:92
TProfile * hResoVSPt_prof_
Definition: Histograms.h:1357
TProfile * hResoVSPt_Ovlap_prof_
Definition: Histograms.h:1362
TProfile * hResoVSPt_Endc_24_prof_
Definition: Histograms.h:1361
TProfile * hResoVSPhi_prof_
Definition: Histograms.h:1367
TProfile * hResoVSPhiMinus_prof_
Definition: Histograms.h:1365
int getXindex(const double &x) const
Definition: Histograms.h:1347
TProfile * hResoVSEta_prof_
Definition: Histograms.h:1363
TProfile * hResoVSPhiPlus_prof_
Definition: Histograms.h:1366
double ** resoVsPtEta_
Definition: Histograms.h:1355
int HFunctionResolution::getXindex ( const double &  x) const
inlineprotected

Definition at line 1347 of file Histograms.h.

References deltaX_, totBinsX_, and xMin_.

Referenced by Fill(), and HFunctionResolutionVarianceCheck::Fill().

1347  {
1348  return int((x-xMin_)/deltaX_*totBinsX_);
1349  }
Definition: DDAxes.h:10
int HFunctionResolution::getYindex ( const double &  y) const
inlineprotected

Definition at line 1350 of file Histograms.h.

References deltaY_, totBinsY_, and yMin_.

Referenced by Fill(), and HFunctionResolutionVarianceCheck::Fill().

1350  {
1351  return int((y-yMin_)/deltaY_*totBinsY_);
1352  }
virtual void HFunctionResolution::Write ( )
inlinevirtual

Implements Histograms.

Reimplemented in HFunctionResolutionVarianceCheck.

Definition at line 1291 of file Histograms.h.

References MultipleCompare::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_.

Referenced by HFunctionResolutionVarianceCheck::Write().

1291  {
1292  if(histoDir_ != 0) histoDir_->cd();
1293 
1294  hReso_->Write();
1295 
1296  for( int i=0; i<totBinsX_; ++i ) {
1297  for( int j=0; j<totBinsY_; ++j ) {
1298  int N = resoCount_[i][j];
1299  // Fill with the mean value
1300  if( N != 0 ) hResoVSPtEta_->SetBinContent( i+1, j+1, resoVsPtEta_[i][j]/N );
1301  else hResoVSPtEta_->SetBinContent( i+1, j+1, 0 );
1302  }
1303  }
1304  hResoVSPtEta_->Write();
1305 
1306  hResoVSPt_prof_->Write();
1307  hResoVSPt_Bar_prof_->Write();
1308  hResoVSPt_Endc_17_prof_->Write();
1309  hResoVSPt_Endc_20_prof_->Write();
1310  hResoVSPt_Endc_24_prof_->Write();
1311  hResoVSPt_Ovlap_prof_->Write();
1312  hResoVSEta_prof_->Write();
1313  //hResoVSTheta_prof_->Write();
1314  hResoVSPhiMinus_prof_->Write();
1315  hResoVSPhiPlus_prof_->Write();
1316  hResoVSPhi_prof_->Write();
1317 
1318  TCanvas canvas(TString(hResoVSPtEta_->GetName())+"_canvas", TString(hResoVSPtEta_->GetTitle())+" canvas", 1000, 800);
1319  canvas.Divide(2);
1320  canvas.cd(1);
1321  hResoVSPtEta_->Draw("lego");
1322  canvas.cd(2);
1323  hResoVSPtEta_->Draw("surf5");
1324  canvas.Write();
1325  hResoVSPtEta_->Write();
1326 
1327  outputFile_->cd();
1328  }
int i
Definition: DBlmapReader.cc:9
TFile * outputFile_
Definition: Histograms.h:107
TProfile * hResoVSPt_Endc_20_prof_
Definition: Histograms.h:1360
TProfile * hResoVSPt_Endc_17_prof_
Definition: Histograms.h:1359
TProfile * hResoVSPt_Bar_prof_
Definition: Histograms.h:1358
TProfile * hResoVSPt_prof_
Definition: Histograms.h:1357
TDirectory * histoDir_
Definition: Histograms.h:108
int j
Definition: DBlmapReader.cc:9
TProfile * hResoVSPt_Ovlap_prof_
Definition: Histograms.h:1362
TProfile * hResoVSPt_Endc_24_prof_
Definition: Histograms.h:1361
TProfile * hResoVSPhi_prof_
Definition: Histograms.h:1367
TProfile * hResoVSPhiMinus_prof_
Definition: Histograms.h:1365
TProfile * hResoVSEta_prof_
Definition: Histograms.h:1363
TProfile * hResoVSPhiPlus_prof_
Definition: Histograms.h:1366
double ** resoVsPtEta_
Definition: Histograms.h:1355

Member Data Documentation

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
int HFunctionResolution::totBinsY_
protected
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().