Adding AliComparisonDraw (Marian)
[u/mrichter/AliRoot.git] / ITS / AliITSPrintRecPoints.C
... / ...
CommitLineData
1void AliITSPrintRecPoints(Int_t outtype=1,TString rfn="galice.root",
2 Int_t mod=-1,Int_t evnt=-1){
3 // Macro to print out the recpoints for all or a specific module
4
5 // Dynamically link some shared libs
6 if (gClassTable->GetID("AliRun") < 0) {
7 gROOT->LoadMacro("loadlibs.C");
8 loadlibs();
9 }
10 else {
11 if(gAlice){
12 delete gAlice->GetRunLoader();
13 delete gAlice;
14 gAlice=0;
15 }
16 }
17
18 gROOT->LoadMacro("$(ALICE_ROOT)/ITS/AliITSstandard.C");
19
20 AliRunLoader *rl = AccessFile(rfn); // Set up to read in Data
21 Int_t retval = rl->LoadHeader();
22 if (retval){
23 cerr<<"AliITSPrintRecPoints.C : LoadHeader returned error"<<endl;
24 return;
25 }
26
27 AliITSLoader* ITSloader = (AliITSLoader*) rl->GetLoader("ITSLoader");
28
29 if(!ITSloader){
30 cerr<<"AliITSPrintRecPoints.C : ITS loader not found"<<endl;
31 return;
32 }
33 cout <<"ITSloader ok"<<endl;
34
35 if(!gGeoManager){
36 gGeoManger = new TGeoManager();
37 gGeoManager->Import("geometry.root","");
38 } // end if
39 ITSloader->LoadHits("read");
40 ITSloader->LoadDigits("read");
41 ITSloader->LoadRecPoints("read");
42 cout << "loaded hits, digits, and RecPoints"<< endl;
43 AliITS *ITS = 0;
44 ITS = (AliITS*)(gAlice->GetDetector("ITS"));
45 if(!ITS){
46 cout << "Error: no ITS found. Aborting"<<endl;
47 return;
48 } // end if !ITS
49 cout <<"ITS="<<ITS<<endl;
50 if(!(ITS->GetDetTypeSim())){
51 cout <<"No AliITSDetTypeSim object found in ITS"<<endl;
52 return;
53 } // end if
54 AliITSgeom *gm=0;
55 gm = ITS->GetITSgeom();
56 if(!gm){
57 cout <<"No AliITSgeom object found in ITS"<<endl;
58 if(!gGeoManager){
59 cout <<"No gGeoManger. Aborting"<<endl;
60 return;
61 }else{
62 ITS->UpdateInternalGeometry(gGeoManager);
63 } // end if
64 } // end if !AliITSgeom
65
66 Int_t evNumber1 = 0;
67 Int_t evNumber2 = gAlice->GetEventsPerRun();
68 if(evnt>=0){
69 evNumber1 = evnt;
70 evNumber2 = evnt+1;
71 } // end if evnt>=0
72 Int_t mod1 = 0;
73 Int_t mod2 = gm->GetIndexMax();
74 if(mod>=0){
75 mod1 = mod;
76 mod2 = mod+1;
77 } // end if mod>=0
78 TClonesArray *rpa;
79 AliITSRecPoint *rp = 0;
80 AliITSDetTypeRec* rec = new AliITSDetTypeRec(ITSloader);
81 //rec->SetITSgeom(gm);
82 rec->SetDefaults();
83
84 Int_t event,m,i,i2;
85 Float_t xyz[3];
86 for(event = evNumber1; event < evNumber2; event++){
87 rl->GetEvent(event);
88 rec->SetTreeAddress();
89 for(m=mod1;m<mod2;m++){
90 rec->ResetRecPoints();
91 TTree *TR = ITSloader->TreeR();
92 TR->GetEvent(m);
93 rpa = rec->RecPoints();
94 i2 = rpa->GetEntriesFast();
95 cout << "Event=" << event << " module=" << m <<
96 " Number of Recpoints=" << i2 <<endl;
97 for(i=0;i<i2;i++){
98 rp = (AliITSRecPoint*)(rpa->At(i));
99 switch(outtype){
100 case 1:
101 rp->GetGlobalXYZ(xyz);
102 cout << i << " lx=" << rp->GetDetLocalX()
103 << " lz=" << rp->GetDetLocalZ()
104 << " x=" << rp->GetX()
105 << " y=" << rp->GetY()<< " z=" << rp->GetZ()
106 <<" gx=" << xyz[0] << " gy="<< xyz[1] <<" gz="<<xyz[2]
107 << endl;
108 break;
109 default:
110 cout << i << " ";
111 rp->Print((ostream*)cout);
112 cout << endl;
113 break;
114 } // end switch
115 } // end for i
116 } // end for m
117 } // end for event
118
119}