FairRoot/PandaRoot
Public Member Functions | Private Member Functions | Private Attributes | List of all members
PndMvdRadDamList Class Reference

#include <PndMvdRadDamList.h>

Public Member Functions

 PndMvdRadDamList (TString fileName)
 
virtual ~PndMvdRadDamList ()
 
Double_t GetWeight (Double_t energy)
 

Private Member Functions

Int_t FindClosestEnergyIndex (Double_t energy)
 
Double_t Interpolate (Double_t energy, Int_t index)
 
 ClassDef (PndMvdRadDamList, 1)
 

Private Attributes

std::vector< std::pair
< Double_t, Double_t > > 
fList
 

Detailed Description

Definition at line 16 of file PndMvdRadDamList.h.

Constructor & Destructor Documentation

PndMvdRadDamList::PndMvdRadDamList ( TString  fileName)

Definition at line 15 of file PndMvdRadDamList.cxx.

References f, fList, i, t, and W.

15  : fList()
16 {
17  Float_t E, W;
18  TFile f(fileName);
19  TTree* t = (TTree*)f.Get("data");
20  t->SetBranchAddress("E_GeV", &E);
21  t->SetBranchAddress("Weight", &W);
22 
23  for (int i = 0; i < t->GetEntriesFast(); i++){
24  t->GetEntry(i);
25  fList.push_back(std::pair<Double_t, Double_t>(E, W));
26  }
27 
28 }
std::vector< std::pair< Double_t, Double_t > > fList
Int_t i
Definition: run_full.C:25
#define W
Definition: createSTT.C:76
TFile * f
Definition: bump_analys.C:12
TTree * t
Definition: bump_analys.C:13
PndMvdRadDamList::~PndMvdRadDamList ( )
virtual

Definition at line 30 of file PndMvdRadDamList.cxx.

31 {
32  // TODO Auto-generated destructor stub
33 }

Member Function Documentation

PndMvdRadDamList::ClassDef ( PndMvdRadDamList  ,
 
)
private
Int_t PndMvdRadDamList::FindClosestEnergyIndex ( Double_t  energy)
private

Definition at line 41 of file PndMvdRadDamList.cxx.

References fList, and i.

Referenced by GetWeight().

42 {
43 
44  for (unsigned int i = 0; i < fList.size(); i++){
45  if (fList[i].first > energy)
46  return i-1;
47  }
48  return -2; //energy is bigger than the latest entry in the list
49 }
std::vector< std::pair< Double_t, Double_t > > fList
Int_t i
Definition: run_full.C:25
Double_t energy
Definition: plot_dirc.C:15
Double_t PndMvdRadDamList::GetWeight ( Double_t  energy)

Definition at line 35 of file PndMvdRadDamList.cxx.

References FindClosestEnergyIndex(), and Interpolate().

36 {
37  Int_t index = FindClosestEnergyIndex(energy);
38  return Interpolate(energy, index);
39 }
Double_t Interpolate(Double_t energy, Int_t index)
Int_t FindClosestEnergyIndex(Double_t energy)
Double_t energy
Definition: plot_dirc.C:15
Double_t PndMvdRadDamList::Interpolate ( Double_t  energy,
Int_t  index 
)
private

Definition at line 51 of file PndMvdRadDamList.cxx.

References fList.

Referenced by GetWeight().

52 {
53  std::pair<Double_t, Double_t> ew1, ew2;
54  if (index == -1){ //energy is lower than the smallest energy entry in the list
55  ew1 = fList[0];
56  ew2 = fList[1];
57  std::cout << "-W- PndMvdRadDamList::Interpolate: Energy is lower than list values" << std::endl;
58  }
59  else if (index < -1){
60  ew1 = fList[fList.size()-2];
61  ew2 = fList[fList.size()-1];
62  std::cout << "-W- PndMvdRadDamList::Interpolate: Energy is bigger than list values" << std::endl;
63  }
64  else{
65  ew1 = fList[index];
66  ew2 = fList[index+1];
67  }
68  return ew1.second + (ew2.second - ew1.second)/(ew2.first - ew1.first) * (energy - ew1.first);
69 }
std::vector< std::pair< Double_t, Double_t > > fList
Double_t energy
Definition: plot_dirc.C:15

Member Data Documentation

std::vector<std::pair<Double_t, Double_t> > PndMvdRadDamList::fList
private

Definition at line 25 of file PndMvdRadDamList.h.

Referenced by FindClosestEnergyIndex(), Interpolate(), and PndMvdRadDamList().


The documentation for this class was generated from the following files: