]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Reding the correction containers in the class
authorpchrist <pchrist@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 5 Aug 2008 10:13:54 +0000 (10:13 +0000)
committerpchrist <pchrist@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 5 Aug 2008 10:13:54 +0000 (10:13 +0000)
PWG2/SPECTRA/AliProtonAnalysis.cxx
PWG2/SPECTRA/AliProtonAnalysis.h

index d2fe151994e18d6bdbfe5f6a881a42a2b18e8069..2cfcc998a5946b872c720716ae905a08855614db 100644 (file)
@@ -70,12 +70,10 @@ AliProtonAnalysis::AliProtonAnalysis() :
   fElectronFunction(0), fMuonFunction(0),
   fPionFunction(0), fKaonFunction(0), fProtonFunction(0),
   fUseTPCOnly(kFALSE), fHistEvents(0), fHistYPtProtons(0), fHistYPtAntiProtons(0),
-  fCorrectionList2D(0), fEfficiencyList1D(0), fCorrectionList1D(0) {
+  fCorrectionListProtons2D(0), fEfficiencyListProtons1D(0), fCorrectionListProtons1D(0),
+  fCorrectionListAntiProtons2D(0), fEfficiencyListAntiProtons1D(0), fCorrectionListAntiProtons1D(0) {
   //Default constructor
   for(Int_t i = 0; i < 5; i++) fPartFrac[i] = 0.0;
-  fCorrectionList2D = new TList(); 
-  fEfficiencyList1D = new TList(); 
-  fCorrectionList1D = new TList();
 }
 
 //____________________________________________________________________//
@@ -107,7 +105,8 @@ AliProtonAnalysis::AliProtonAnalysis(Int_t nbinsY, Float_t fLowY, Float_t fHighY
   fElectronFunction(0), fMuonFunction(0),
   fPionFunction(0), fKaonFunction(0), fProtonFunction(0),
   fUseTPCOnly(kFALSE), fHistEvents(0), fHistYPtProtons(0), fHistYPtAntiProtons(0),
-  fCorrectionList2D(0), fEfficiencyList1D(0), fCorrectionList1D(0) {
+  fCorrectionListProtons2D(0), fEfficiencyListProtons1D(0), fCorrectionListProtons1D(0),
+  fCorrectionListAntiProtons2D(0), fEfficiencyListAntiProtons1D(0), fCorrectionListAntiProtons1D(0) {
   //Default constructor
 
   fHistEvents = new TH1I("fHistEvents","Analyzed events",1,0,1);
@@ -131,9 +130,6 @@ AliProtonAnalysis::~AliProtonAnalysis() {
   if(fHistEvents) delete fHistEvents;
   if(fHistYPtProtons) delete fHistYPtProtons;
   if(fHistYPtAntiProtons) delete fHistYPtAntiProtons;
-  if(fCorrectionList2D) delete fCorrectionList2D;
-  if(fEfficiencyList1D) delete fEfficiencyList1D;
-  if(fCorrectionList1D) delete fCorrectionList1D;
   if(fGlobalQAList) delete fGlobalQAList;
   if(fQA2DList) delete fQA2DList;
   if(fQAPrimaryProtonsAcceptedList) delete fQAPrimaryProtonsAcceptedList;
@@ -144,6 +140,12 @@ AliProtonAnalysis::~AliProtonAnalysis() {
   if(fQAPrimaryAntiProtonsRejectedList) delete fQAPrimaryAntiProtonsRejectedList;
   if(fQASecondaryAntiProtonsAcceptedList) delete fQASecondaryAntiProtonsAcceptedList;
   if(fQASecondaryAntiProtonsRejectedList) delete fQASecondaryAntiProtonsRejectedList; 
+  if(fCorrectionListProtons2D) delete fCorrectionListProtons2D;
+  if(fEfficiencyListProtons1D) delete fEfficiencyListProtons1D;
+  if(fCorrectionListProtons1D) delete fCorrectionListProtons1D;
+  if(fCorrectionListAntiProtons2D) delete fCorrectionListAntiProtons2D;
+  if(fEfficiencyListAntiProtons1D) delete fEfficiencyListAntiProtons1D;
+  if(fCorrectionListAntiProtons1D) delete fCorrectionListAntiProtons1D;
 }
 
 //____________________________________________________________________//
@@ -1221,19 +1223,25 @@ Bool_t AliProtonAnalysis::ReadCorrectionContainer(const char* filename) {
     status = kFALSE;
   }
 
-  AliCFContainer *corrfwContainer = (AliCFContainer*) (file->Get("container"));
-  if(!corrfwContainer) {
-    cout<<"CORRFW container not found!"<<endl;
+  //________________________________________//
+  //Protons
+  fCorrectionListProtons2D = new TList(); 
+  fEfficiencyListProtons1D = new TList(); 
+  fCorrectionListProtons1D = new TList();
+  
+  AliCFContainer *corrfwContainerProtons = (AliCFContainer*) (file->Get("containerProtons"));
+  if(!corrfwContainerProtons) {
+    cout<<"CORRFW container for protons not found!"<<endl;
     status = kFALSE;
   }
   
-  Int_t nSteps = corrfwContainer->GetNStep();
+  Int_t nSteps = corrfwContainerProtons->GetNStep();
   TH2D *gYPt[4];
   //currently the GRID is formed by the y-pT parameters
   //Add Vz as a next step
   Int_t iRap = 0, iPt = 1;
   for(Int_t iStep = 0; iStep < nSteps; iStep++) {
-    gYPt[iStep] = corrfwContainer->ShowProjection(iRap,iPt,iStep);
+    gYPt[iStep] = corrfwContainerProtons->ShowProjection(iRap,iPt,iStep);
     //fCorrectionList2D->Add(gYPt[iStep]);
   }
 
@@ -1245,18 +1253,71 @@ Bool_t AliProtonAnalysis::ReadCorrectionContainer(const char* filename) {
 
   //Get the 2D efficiency maps
   for(Int_t iStep = 1; iStep < nSteps; iStep++) { 
-    gTitle = "Efficiency_Step0_Step"; gTitle += iStep; 
+    gTitle = "EfficiencyProtons_Step0_Step"; gTitle += iStep; 
+    efficiency[iStep] = new AliCFEffGrid(gTitle.Data(),
+                                        gTitle.Data(),*corrfwContainerProtons);
+    efficiency[iStep]->CalculateEfficiency(iStep,0); //eff= step[i]/step0
+    fCorrectionListProtons2D->Add(efficiency[iStep]);  
+  }
+  //Get the projection of the efficiency maps
+  for(Int_t iParameter = 0; iParameter < 2; iParameter++) { 
+    for(Int_t iStep = 1; iStep < nSteps; iStep++) { 
+      gEfficiency[iParameter][iStep-1] = efficiency[iStep]->Project(iParameter);
+      fEfficiencyListProtons1D->Add(gEfficiency[iParameter][iStep-1]);  
+      gTitle = "ProtonsCorrection_Parameter"; gTitle += iParameter+1;
+      gTitle += "_Step0_Step"; gTitle += iStep; 
+      gCorrection[iParameter][iStep-1] = new TH1D(gTitle.Data(),
+                                                  gTitle.Data(),
+                                                  gEfficiency[iParameter][iStep-1]->GetNbinsX(),
+                                                  gEfficiency[iParameter][iStep-1]->GetXaxis()->GetXmin(),
+                                                  gEfficiency[iParameter][iStep-1]->GetXaxis()->GetXmax());
+      //initialisation of the correction
+      for(Int_t iBin = 1; iBin <= gEfficiency[iParameter][iStep-1]->GetNbinsX(); iBin++)
+       gCorrection[iParameter][iStep-1]->SetBinContent(iBin,1.0);
+    }//step loop
+  }//parameter loop
+  //Calculate the 1D correction parameters as a function of y and pT
+  for(Int_t iParameter = 0; iParameter < 2; iParameter++) { 
+    for(Int_t iStep = 1; iStep < nSteps; iStep++) { 
+      gCorrection[iParameter][iStep-1]->Divide(gEfficiency[iParameter][iStep-1]);
+      fCorrectionListProtons1D->Add(gCorrection[iParameter][iStep-1]);  
+    }
+  }
+
+  //________________________________________//
+  //AntiProtons
+  fCorrectionListAntiProtons2D = new TList(); 
+  fEfficiencyListAntiProtons1D = new TList(); 
+  fCorrectionListAntiProtons1D = new TList();
+  
+  AliCFContainer *corrfwContainerAntiProtons = (AliCFContainer*) (file->Get("containerAntiProtons"));
+  if(!corrfwContainerAntiProtons) {
+    cout<<"CORRFW container for antiprotons not found!"<<endl;
+    status = kFALSE;
+  }
+  
+  nSteps = corrfwContainerAntiProtons->GetNStep();
+  //currently the GRID is formed by the y-pT parameters
+  //Add Vz as a next step
+  iRap = 0; iPt = 1;
+  for(Int_t iStep = 0; iStep < nSteps; iStep++) {
+    gYPt[iStep] = corrfwContainerAntiProtons->ShowProjection(iRap,iPt,iStep);
+  }
+
+  //Get the 2D efficiency maps
+  for(Int_t iStep = 1; iStep < nSteps; iStep++) { 
+    gTitle = "EfficiencyAntiProtons_Step0_Step"; gTitle += iStep; 
     efficiency[iStep] = new AliCFEffGrid(gTitle.Data(),
-                                        gTitle.Data(),*corrfwContainer);
+                                        gTitle.Data(),*corrfwContainerAntiProtons);
     efficiency[iStep]->CalculateEfficiency(iStep,0); //eff= step[i]/step0
-    fCorrectionList2D->Add(efficiency[iStep]);  
+    fCorrectionListAntiProtons2D->Add(efficiency[iStep]);  
   }
   //Get the projection of the efficiency maps
   for(Int_t iParameter = 0; iParameter < 2; iParameter++) { 
     for(Int_t iStep = 1; iStep < nSteps; iStep++) { 
       gEfficiency[iParameter][iStep-1] = efficiency[iStep]->Project(iParameter);
-      fEfficiencyList1D->Add(gEfficiency[iParameter][iStep-1]);  
-      gTitle = "Correction_Parameter"; gTitle += iParameter+1;
+      fEfficiencyListProtons1D->Add(gEfficiency[iParameter][iStep-1]);  
+      gTitle = "AntiProtonsCorrection_Parameter"; gTitle += iParameter+1;
       gTitle += "_Step0_Step"; gTitle += iStep; 
       gCorrection[iParameter][iStep-1] = new TH1D(gTitle.Data(),
                                                   gTitle.Data(),
@@ -1272,7 +1333,7 @@ Bool_t AliProtonAnalysis::ReadCorrectionContainer(const char* filename) {
   for(Int_t iParameter = 0; iParameter < 2; iParameter++) { 
     for(Int_t iStep = 1; iStep < nSteps; iStep++) { 
       gCorrection[iParameter][iStep-1]->Divide(gEfficiency[iParameter][iStep-1]);
-      fCorrectionList1D->Add(gCorrection[iParameter][iStep-1]);  
+      fCorrectionListAntiProtons1D->Add(gCorrection[iParameter][iStep-1]);  
     }
   }
 
index 72e6393e9540497143c2c301ba747a8c235ca42e..1874804f6038c3fee0fb40b92452e42c0fb78a64 100644 (file)
@@ -137,9 +137,12 @@ class AliProtonAnalysis : public TObject {
 
   //interface to the correction framework
   Bool_t ReadCorrectionContainer(const char* filename);
-  TList *GetCorrectionList2D() {return fCorrectionList2D;} 
-  TList *GetEfficiencyList1D() {return fEfficiencyList1D;} 
-  TList *GetCorrectionList1D() {return fCorrectionList1D;} 
+  TList *GetCorrectionListProtons2D() {return fCorrectionListProtons2D;} 
+  TList *GetEfficiencyListProtons1D() {return fEfficiencyListProtons1D;} 
+  TList *GetCorrectionListProtons1D() {return fCorrectionListProtons1D;} 
+  TList *GetCorrectionListAntiProtons2D() {return fCorrectionListAntiProtons2D;} 
+  TList *GetEfficiencyListAntiProtons1D() {return fEfficiencyListAntiProtons1D;} 
+  TList *GetCorrectionListAntiProtons1D() {return fCorrectionListAntiProtons1D;} 
 
  private:
   AliProtonAnalysis(const AliProtonAnalysis&); // Not implemented
@@ -181,114 +184,6 @@ class AliProtonAnalysis : public TObject {
   TList *fQAPrimaryAntiProtonsRejectedList; //list of the QA histos for rejected primary antiprotons
   TList *fQASecondaryAntiProtonsAcceptedList; //list of the QA histos for accepted secondary antiprotons
   TList *fQASecondaryAntiProtonsRejectedList; //list of the QA histos for rejected secondary antiprotons
-  //primary protons
-  /*TH1F *fPrimaryProtonsTPCClustersReject;        //QA histogram for the primaries rejected by the TPC cluster cut
-  TH1F *fPrimaryProtonsTPCClustersPass;          //QA histogram for the primaries accepted by the TPC cluster cut
-  TH1F *fPrimaryProtonsChi2PerClusterTPCReject;  //QA histogram for the primaries rejected by the chi2 per TPC cluster cut 
-  TH1F *fPrimaryProtonsChi2PerClusterTPCPass;    //QA histogram for the primaries accepted by the chi2 per TPC cluster cut
-  TH1F *fPrimaryProtonsExtCov11Reject;           //QA histogram for the primaries rejected by the sigma of the local Y cut
-  TH1F *fPrimaryProtonsExtCov11Pass;             //QA histogram for the primaries accepted by the sigma of the local Y cut
-  TH1F *fPrimaryProtonsExtCov22Reject;           //QA histogram for the primaries rejected by the sigma of the local Z cut 
-  TH1F *fPrimaryProtonsExtCov22Pass;             //QA histogram for the primaries accepted by the sigma of the local Z cut 
-  TH1F *fPrimaryProtonsExtCov33Reject;           //QA histogram for the primaries rejected by the sigma of the sin(phi) cut 
-  TH1F *fPrimaryProtonsExtCov33Pass;             //QA histogram for the primaries accepted by the sigma of the sin(phi) cut 
-  TH1F *fPrimaryProtonsExtCov44Reject;           //QA histogram for the primaries rejected by the sigma of the tan(lambda) cut 
-  TH1F *fPrimaryProtonsExtCov44Pass;             //QA histogram for the primaries accepted by the sigma of the tan(lambda) cut 
-  TH1F *fPrimaryProtonsExtCov55Reject;           //QA histogram for the primaries rejected by the the sigma of 1/pT cut 
-  TH1F *fPrimaryProtonsExtCov55Pass;             //QA histogram for the primaries accepted by the sigma of the 1/pT cut 
-  TH1F *fPrimaryProtonsSigmaToVertexReject;      //QA histogram for the primaries rejected by the sigma to vertex cut 
-  TH1F *fPrimaryProtonsSigmaToVertexPass;        //QA histogram for the primaries accepted by the sigma to vertex cut 
-  TH1F *fPrimaryProtonsSigmaToVertexTPCReject;   //QA histogram for the primaries rejected by the sigma to vertex (TPC) cut 
-  TH1F *fPrimaryProtonsSigmaToVertexTPCPass;     //QA histogram for the primaries accepted by the sigma to vertex (TPC) cut 
-  TH1F *fPrimaryProtonsITSRefitReject;           //QA histogram for the primaries rejected by the ITS refit cut 
-  TH1F *fPrimaryProtonsITSRefitPass;             //QA histogram for the primaries accepted by the ITS refit cut
-  TH1F *fPrimaryProtonsTPCRefitReject;           //QA histogram for the primaries rejected by the TPC refit cut
-  TH1F *fPrimaryProtonsTPCRefitPass;             //QA histogram for the primaries accepted by the TPC refit cut
-  TH1F *fPrimaryProtonsESDpidReject;             //QA histogram for the primaries rejected by the ESD pid cut
-  TH1F *fPrimaryProtonsESDpidPass;               //QA histogram for the primaries accepted by the ESD pid cut
-  TH1F *fPrimaryProtonsTPCpidReject;             //QA histogram for the primaries rejected by the TPC pid cut
-  TH1F *fPrimaryProtonsTPCpidPass;               //QA histogram for the primaries accepted by the TPC pid cut
-  //secondary protons
-  TH1F *fSecondaryProtonsTPCClustersReject;        //QA histogram for the secondaries rejected by the TPC cluster cut
-  TH1F *fSecondaryProtonsTPCClustersPass;          //QA histogram for the secondaries accepted by the TPC cluster cut
-  TH1F *fSecondaryProtonsChi2PerClusterTPCReject;  //QA histogram for the secondaries rejected by the chi2 per TPC cluster cut 
-  TH1F *fSecondaryProtonsChi2PerClusterTPCPass;    //QA histogram for the secondaries accepted by the chi2 per TPC cluster cut
-  TH1F *fSecondaryProtonsExtCov11Reject;           //QA histogram for the secondaries rejected by the sigma of the local Y cut
-  TH1F *fSecondaryProtonsExtCov11Pass;             //QA histogram for the secondaries accepted by the sigma of the local Y cut
-  TH1F *fSecondaryProtonsExtCov22Reject;           //QA histogram for the secondaries rejected by the sigma of the local Z cut 
-  TH1F *fSecondaryProtonsExtCov22Pass;             //QA histogram for the secondaries accepted by the sigma of the local Z cut 
-  TH1F *fSecondaryProtonsExtCov33Reject;           //QA histogram for the secondaries rejected by the sigma of the sin(phi) cut 
-  TH1F *fSecondaryProtonsExtCov33Pass;             //QA histogram for the secondaries accepted by the sigma of the sin(phi) cut 
-  TH1F *fSecondaryProtonsExtCov44Reject;           //QA histogram for the secondaries rejected by the sigma of the tan(lambda) cut 
-  TH1F *fSecondaryProtonsExtCov44Pass;             //QA histogram for the secondaries accepted by the sigma of the tan(lambda) cut 
-  TH1F *fSecondaryProtonsExtCov55Reject;           //QA histogram for the secondaries rejected by the the sigma of 1/pT cut 
-  TH1F *fSecondaryProtonsExtCov55Pass;             //QA histogram for the secondaries accepted by the sigma of the 1/pT cut 
-  TH1F *fSecondaryProtonsSigmaToVertexReject;      //QA histogram for the secondaries rejected by the sigma to vertex cut 
-  TH1F *fSecondaryProtonsSigmaToVertexPass;        //QA histogram for the secondaries accepted by the sigma to vertex cut 
-  TH1F *fSecondaryProtonsSigmaToVertexTPCReject;   //QA histogram for the secondaries rejected by the sigma to vertex (TPC) cut 
-  TH1F *fSecondaryProtonsSigmaToVertexTPCPass;     //QA histogram for the secondaries accepted by the sigma to vertex (TPC) cut 
-  TH1F *fSecondaryProtonsITSRefitReject;           //QA histogram for the secondaries rejected by the ITS refit cut 
-  TH1F *fSecondaryProtonsITSRefitPass;             //QA histogram for the secondaries accepted by the ITS refit cut
-  TH1F *fSecondaryProtonsTPCRefitReject;           //QA histogram for the secondaries rejected by the TPC refit cut
-  TH1F *fSecondaryProtonsTPCRefitPass;             //QA histogram for the secondaries accepted by the TPC refit cut
-  TH1F *fSecondaryProtonsESDpidReject;             //QA histogram for the secondaries rejected by the ESD pid cut
-  TH1F *fSecondaryProtonsESDpidPass;               //QA histogram for the secondaries accepted by the ESD pid cut
-  TH1F *fSecondaryProtonsTPCpidReject;             //QA histogram for the secondaries rejected by the TPC pid cut
-  TH1F *fSecondaryProtonsTPCpidPass;               //QA histogram for the secondaries accepted by the TPC pid cut
-  //primary antiprotons
-  TH1F *fPrimaryAntiProtonsTPCClustersReject;        //QA histogram for the primaries rejected by the TPC cluster cut
-  TH1F *fPrimaryAntiProtonsTPCClustersPass;          //QA histogram for the primaries accepted by the TPC cluster cut
-  TH1F *fPrimaryAntiProtonsChi2PerClusterTPCReject;  //QA histogram for the primaries rejected by the chi2 per TPC cluster cut 
-  TH1F *fPrimaryAntiProtonsChi2PerClusterTPCPass;    //QA histogram for the primaries accepted by the chi2 per TPC cluster cut
-  TH1F *fPrimaryAntiProtonsExtCov11Reject;           //QA histogram for the primaries rejected by the sigma of the local Y cut
-  TH1F *fPrimaryAntiProtonsExtCov11Pass;             //QA histogram for the primaries accepted by the sigma of the local Y cut
-  TH1F *fPrimaryAntiProtonsExtCov22Reject;           //QA histogram for the primaries rejected by the sigma of the local Z cut 
-  TH1F *fPrimaryAntiProtonsExtCov22Pass;             //QA histogram for the primaries accepted by the sigma of the local Z cut 
-  TH1F *fPrimaryAntiProtonsExtCov33Reject;           //QA histogram for the primaries rejected by the sigma of the sin(phi) cut 
-  TH1F *fPrimaryAntiProtonsExtCov33Pass;             //QA histogram for the primaries accepted by the sigma of the sin(phi) cut 
-  TH1F *fPrimaryAntiProtonsExtCov44Reject;           //QA histogram for the primaries rejected by the sigma of the tan(lambda) cut 
-  TH1F *fPrimaryAntiProtonsExtCov44Pass;             //QA histogram for the primaries accepted by the sigma of the tan(lambda) cut 
-  TH1F *fPrimaryAntiProtonsExtCov55Reject;           //QA histogram for the primaries rejected by the the sigma of 1/pT cut 
-  TH1F *fPrimaryAntiProtonsExtCov55Pass;             //QA histogram for the primaries accepted by the sigma of the 1/pT cut 
-  TH1F *fPrimaryAntiProtonsSigmaToVertexReject;      //QA histogram for the primaries rejected by the sigma to vertex cut 
-  TH1F *fPrimaryAntiProtonsSigmaToVertexPass;        //QA histogram for the primaries accepted by the sigma to vertex cut 
-  TH1F *fPrimaryAntiProtonsSigmaToVertexTPCReject;   //QA histogram for the primaries rejected by the sigma to vertex (TPC) cut 
-  TH1F *fPrimaryAntiProtonsSigmaToVertexTPCPass;     //QA histogram for the primaries accepted by the sigma to vertex (TPC) cut 
-  TH1F *fPrimaryAntiProtonsITSRefitReject;           //QA histogram for the primaries rejected by the ITS refit cut 
-  TH1F *fPrimaryAntiProtonsITSRefitPass;             //QA histogram for the primaries accepted by the ITS refit cut
-  TH1F *fPrimaryAntiProtonsTPCRefitReject;           //QA histogram for the primaries rejected by the TPC refit cut
-  TH1F *fPrimaryAntiProtonsTPCRefitPass;             //QA histogram for the primaries accepted by the TPC refit cut
-  TH1F *fPrimaryAntiProtonsESDpidReject;             //QA histogram for the primaries rejected by the ESD pid cut
-  TH1F *fPrimaryAntiProtonsESDpidPass;               //QA histogram for the primaries accepted by the ESD pid cut
-  TH1F *fPrimaryAntiProtonsTPCpidReject;             //QA histogram for the primaries rejected by the TPC pid cut
-  TH1F *fPrimaryAntiProtonsTPCpidPass;               //QA histogram for the primaries accepted by the TPC pid cut
-  //secondary antiprotons
-  TH1F *fSecondaryAntiProtonsTPCClustersReject;        //QA histogram for the secondaries rejected by the TPC cluster cut
-  TH1F *fSecondaryAntiProtonsTPCClustersPass;          //QA histogram for the secondaries accepted by the TPC cluster cut
-  TH1F *fSecondaryAntiProtonsChi2PerClusterTPCReject;  //QA histogram for the secondaries rejected by the chi2 per TPC cluster cut 
-  TH1F *fSecondaryAntiProtonsChi2PerClusterTPCPass;    //QA histogram for the secondaries accepted by the chi2 per TPC cluster cut
-  TH1F *fSecondaryAntiProtonsExtCov11Reject;           //QA histogram for the secondaries rejected by the sigma of the local Y cut
-  TH1F *fSecondaryAntiProtonsExtCov11Pass;             //QA histogram for the secondaries accepted by the sigma of the local Y cut
-  TH1F *fSecondaryAntiProtonsExtCov22Reject;           //QA histogram for the secondaries rejected by the sigma of the local Z cut 
-  TH1F *fSecondaryAntiProtonsExtCov22Pass;             //QA histogram for the secondaries accepted by the sigma of the local Z cut 
-  TH1F *fSecondaryAntiProtonsExtCov33Reject;           //QA histogram for the secondaries rejected by the sigma of the sin(phi) cut 
-  TH1F *fSecondaryAntiProtonsExtCov33Pass;             //QA histogram for the secondaries accepted by the sigma of the sin(phi) cut 
-  TH1F *fSecondaryAntiProtonsExtCov44Reject;           //QA histogram for the secondaries rejected by the sigma of the tan(lambda) cut 
-  TH1F *fSecondaryAntiProtonsExtCov44Pass;             //QA histogram for the secondaries accepted by the sigma of the tan(lambda) cut 
-  TH1F *fSecondaryAntiProtonsExtCov55Reject;           //QA histogram for the secondaries rejected by the the sigma of 1/pT cut 
-  TH1F *fSecondaryAntiProtonsExtCov55Pass;             //QA histogram for the secondaries accepted by the sigma of the 1/pT cut 
-  TH1F *fSecondaryAntiProtonsSigmaToVertexReject;      //QA histogram for the secondaries rejected by the sigma to vertex cut 
-  TH1F *fSecondaryAntiProtonsSigmaToVertexPass;        //QA histogram for the secondaries accepted by the sigma to vertex cut 
-  TH1F *fSecondaryAntiProtonsSigmaToVertexTPCReject;   //QA histogram for the secondaries rejected by the sigma to vertex (TPC) cut 
-  TH1F *fSecondaryAntiProtonsSigmaToVertexTPCPass;     //QA histogram for the secondaries accepted by the sigma to vertex (TPC) cut 
-  TH1F *fSecondaryAntiProtonsITSRefitReject;           //QA histogram for the secondaries rejected by the ITS refit cut 
-  TH1F *fSecondaryAntiProtonsITSRefitPass;             //QA histogram for the secondaries accepted by the ITS refit cut
-  TH1F *fSecondaryAntiProtonsTPCRefitReject;           //QA histogram for the secondaries rejected by the TPC refit cut
-  TH1F *fSecondaryAntiProtonsTPCRefitPass;             //QA histogram for the secondaries accepted by the TPC refit cut
-  TH1F *fSecondaryAntiProtonsESDpidReject;             //QA histogram for the secondaries rejected by the ESD pid cut
-  TH1F *fSecondaryAntiProtonsESDpidPass;               //QA histogram for the secondaries accepted by the ESD pid cut
-  TH1F *fSecondaryAntiProtonsTPCpidReject;             //QA histogram for the secondaries rejected by the TPC pid cut
-  TH1F *fSecondaryAntiProtonsTPCpidPass;*/               //QA histogram for the secondaries accepted by the TPC pid cut
 
   //pid
   Bool_t fFunctionProbabilityFlag; //flag: kTRUE if functions used
@@ -307,9 +202,12 @@ class AliProtonAnalysis : public TObject {
   TH2F *fHistYPtAntiProtons; // Y-Pt of Antiprotons
 
   //Corrections
-  TList *fCorrectionList2D; //list for the 2d corrections 
-  TList *fEfficiencyList1D; //list for the 1d efficiencies
-  TList *fCorrectionList1D; //list for the 1d corrections 
+  TList *fCorrectionListProtons2D; //list for the 2d corrections 
+  TList *fEfficiencyListProtons1D; //list for the 1d efficiencies
+  TList *fCorrectionListProtons1D; //list for the 1d corrections 
+  TList *fCorrectionListAntiProtons2D; //list for the 2d corrections 
+  TList *fEfficiencyListAntiProtons1D; //list for the 1d efficiencies
+  TList *fCorrectionListAntiProtons1D; //list for the 1d corrections 
   
   ClassDef(AliProtonAnalysis,0);
 };