]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSVertexerZTest.C
Put index into track's name instead of its label.
[u/mrichter/AliRoot.git] / ITS / AliITSVertexerZTest.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include <TFile.h>
3 #include <TGeoManager.h>
4 #include <TH1F.h>
5 #include <TH2F.h>
6 #include <TTree.h>
7 #include <Riostream.h>
8 #include <AliRun.h>
9 #include <AliHeader.h>
10 #include <AliGenEventHeader.h>
11 #include <AliGeomManager.h>
12 #include <AliITSVertexerZ.h>
13 #include <AliRunLoader.h>
14 #include <AliITSLoader.h>
15
16 #endif
17
18 void AliITSVertexerZTest(Float_t delphi=0.05,Float_t window=3.,Float_t initx=0., Float_t inity=0.,TString FileIn="galice.root"){
19   // delphi ---> azimuthal range to accept tracklets
20   // window ---> window in Z around the peak of tracklets proj. in mm
21   Int_t kDebug = 50;
22   TH2F *diff2 = new TH2F("diff2","Zfound vs Ztrue",100,-6,6,100,-6,6); 
23   TH1F *diff1 = new TH1F("diff1","Zfound - Ztrue(#mu m)",100,-500,500);
24   delete gAlice;
25   gAlice = 0;
26   AliRunLoader *rl = AliRunLoader::Open(FileIn.Data());
27   Int_t retval = rl->LoadgAlice();
28   if(retval){
29     cerr<<"AliITSVertexerZTest.C: AliRun object not found"<<endl;
30     return;
31   }
32   retval = rl->LoadHeader();
33   if (retval){
34     cerr<<"AliITSVertexerZTest.C : LoadHeader returned error"<<endl;
35     return;
36   }
37   retval = rl->LoadKinematics();
38   if (retval){
39     cerr<<"AliITSVertexerZTest.C : LoadKinematics returned error"<<endl;
40     return;
41   }
42   
43   AliGeomManager::LoadGeometry("geometry.root");
44
45   AliITSVertexerZ *dovert = new AliITSVertexerZ("default",initx,inity);
46   //dovert->SetDebug(0);
47   //  dovert->SetDiffPhiMax(delphi);
48   //  dovert->SetWindow(window);
49   dovert->PrintStatus();
50   Int_t sigmazero=0;
51   AliESDVertex *vert = 0;
52   for(Int_t i=0; i<rl->TreeE()->GetEntries(); i++){
53     rl->GetEvent(i);
54     // The true Z coord. is fetched for comparison
55     AliHeader *header = rl->GetHeader();
56     AliGenEventHeader* genEventHeader = header->GenEventHeader();
57     TArrayF primaryVertex(3);
58     genEventHeader->PrimaryVertex(primaryVertex);
59     vert = dovert->FindVertexForCurrentEvent(i);
60     if(kDebug>0){
61       // Prints the results
62       cout <<"========================================================\n";
63       cout << "Event number: "<<i<<")  Z Vertex:"<<endl;
64       if(vert){
65         cout<<"FOUND: "<<vert->GetZv()<<"; ";
66         cout<<vert->GetZRes()<<"; "<<vert->GetNContributors()<<endl;
67         cout <<" True Z position "<<primaryVertex[2]<<", diff= ";
68         cout<<(primaryVertex[2]-vert->GetZv())*10000.<<endl;
69       } else {
70         cout<<"NOT FOUND"<<endl;
71       }
72     }
73     if(vert){
74       Float_t found = vert->GetZv();
75       diff2->Fill(primaryVertex[2],found);
76       found = 10000.*(found-primaryVertex[2]);
77       if(vert->GetZRes()!=0){
78         diff1->Fill(found);
79       } else {
80         sigmazero++;
81       }
82       dovert->WriteCurrentVertex();
83     }
84   }
85   if(kDebug>0){
86     cout<<"Only one tracklet (sigma = 0) "<<sigmazero<<endl;
87   }
88   
89 }