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