]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - ITS/AliITSComparisonV2.C
Bayesian PID: new parametrization and code update (E. Biolcati - F.Prino)
[u/mrichter/AliRoot.git] / ITS / AliITSComparisonV2.C
index e142816392b750b95978776ab9cf03726a217c17..34aa7df27772fca536a10602164ba1e57d997bcf 100644 (file)
   #include "AliTrackReference.h"
   #include "AliRunLoader.h"
   #include "AliRun.h"
-  #include "AliESD.h"
+  #include "AliESDEvent.h"
+  #include "AliESDtrack.h"
 
-  #include "AliITS.h"
-  #include "AliITSgeom.h"
   #include "AliITSRecPoint.h"
   #include "AliITSLoader.h"
 #endif
@@ -48,7 +47,7 @@ static Int_t allselected=0;
 static Int_t allfound=0;
 
 Int_t AliITSComparisonV2
-(Float_t ptcutl=0.2, Float_t ptcuth=10., const Char_t *dir=".", Float_t ratio=0.0) {
+(Float_t ptcutl=0.2, Float_t ptcuth=10., const Char_t *dir=".") {
    gBenchmark->Start("AliITSComparisonV2");
 
    ::Info("AliITSComparisonV2.C","Doing comparison...");
@@ -95,7 +94,7 @@ Int_t AliITSComparisonV2
 
    TH1F *he=(TH1F*)gROOT->FindObject("he");
    if (!he) 
-      he =new TH1F("he","dE/dX for pions with 0.4<p<0.5 GeV/c",50,0.,100.);
+      he =new TH1F("he","dE/dX for pions with 0.4<p<0.5 GeV/c",50,0.,200.);
 
    TH2F *hep=(TH2F*)gROOT->FindObject("hep");
    if (!hep) hep=new TH2F("hep","dE/dX vs momentum",50,0.,2.,50,0.,400.);
@@ -144,13 +143,13 @@ Int_t AliITSComparisonV2
          return 4;
       }
    }
-   AliESD* event = new AliESD;
+   AliESDEvent* event = new AliESDEvent();
    TTree* esdTree = (TTree*) ef->Get("esdTree");
    if (!esdTree) {
       ::Error("AliITSComparison.C", "no ESD tree found");
       return 6;
    }
-   esdTree->SetBranchAddress("ESD", &event);
+   event->ReadFromTree(esdTree);
 
 
    //******* Loop over events *********
@@ -208,8 +207,6 @@ Int_t AliITSComparisonV2
            numb[nmult]=cnt; nmult++;        
         }
 
-       if (esd->GetITSFakeRatio()<ratio) tlab = TMath::Abs(tlab);
-
         if (lab==tlab) hfound->Fill(ptg);
         else {
           fake[nfake++]=lab;
@@ -271,6 +268,7 @@ Int_t AliITSComparisonV2
    } //***** End of the loop over events
 
    delete event;
+   delete esdTree;
    ef->Close();
    
    delete itsTree;
@@ -366,7 +364,7 @@ Int_t AliITSComparisonV2
 
 Int_t GoodTracksITS(const Char_t *dir) {
    if (gAlice) { 
-       delete gAlice->GetRunLoader();
+       delete AliRunLoader::Instance();
        delete gAlice;//if everything was OK here it is already NULL
        gAlice = 0x0;
    }
@@ -384,19 +382,6 @@ Int_t GoodTracksITS(const Char_t *dir) {
    rl->LoadHeader();
    rl->LoadKinematics();
 
-   AliITS *ITS=(AliITS*)rl->GetAliRun()->GetDetector("ITS");
-   if (!ITS) {
-      ::Error("GoodTracksITS","Can't get the ITS !");
-      delete rl;
-      return 2;
-   }
-   AliITSgeom *geom=ITS->GetITSgeom();
-   if (!geom) {
-      ::Error("GoodTracksITS","Can't get the ITS geometry !"); 
-      delete rl;
-      return 3;
-   }
-
    AliITSLoader* itsl = (AliITSLoader*)rl->GetLoader("ITSLoader");
    if (itsl == 0x0) {
        ::Error("GoodTracksITS","Can not find the ITSLoader");
@@ -468,14 +453,16 @@ Int_t GoodTracksITS(const Char_t *dir) {
      for (k=0; k<entr; k++) {
          cTree->GetEvent(k);
          Int_t ncl=clusters->GetEntriesFast(); if (ncl==0) continue;
-         Int_t lay,lad,det;  geom->GetModuleId(k,lay,lad,det);
-         if (lay<1 || lay>6) {
-           ::Error("GoodTracksITS","Wrong layer !"); 
-            delete rl;
-            return 10;
-         }
          while (ncl--) {
             AliITSRecPoint *pnt=(AliITSRecPoint*)clusters->UncheckedAt(ncl);
+
+            Int_t lay=pnt->GetLayer();
+            if (lay<0 || lay>5) {
+               ::Error("GoodTracksITS","Wrong layer !");
+               delete rl;
+               return 10;
+            }
+
             Int_t l0=pnt->GetLabel(0);
               if (l0>=np) {
 //              cerr<<"Wrong label: "<<l0<<endl;
@@ -491,7 +478,7 @@ Int_t GoodTracksITS(const Char_t *dir) {
 //              cerr<<"Wrong label: "<<l2<<endl;
                 continue;
               }
-            Int_t mask=1<<(lay-1);
+            Int_t mask=1<<lay;
             if (l0>=0) good[l0]|=mask; 
             if (l1>=0) good[l1]|=mask; 
             if (l2>=0) good[l2]|=mask;