]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/TOFquickanal.C
Updated version of the HMPID DA. To be checked on real data by the experts and then...
[u/mrichter/AliRoot.git] / TOF / TOFquickanal.C
1 Int_t TOFquickanal(Int_t eventNumber = 0)
2 {
3   /////////////////////////////////////////////////////////////////////////
4   //   This macro is a small example of a ROOT macro
5   //   illustrating how to read the output of GALICE
6   //   and fill some histograms concerning the TOF Hit Tree.
7   //
8   //     Root > .L TOFquickanal.C   //this loads the macro in memory
9   //     Root > TOFquickanal();     //by default process first event
10   //     Root > TOFquickanal(2);    //process third event
11   //Begin_Html
12   /*
13     <img src="picts/TOFquickanal.gif">
14   */
15   //End_Html
16   //
17   // Author: F. Pierella , Bologna University 12-04-2001
18   // Updated to the new I/O by: A. De Caro, C. Zampolli
19   /////////////////////////////////////////////////////////////////////////
20   
21   // Dynamically link some shared libs
22   if (gClassTable->GetID("AliRun") < 0) {
23     gROOT->LoadMacro("loadlibs.C");
24     loadlibs();
25   }
26
27   Int_t rc = 0;
28   
29   AliRunLoader *rl =AliRunLoader::Open("galice.root",AliConfig::GetDefaultEventFolderName(),"update");
30   if (!rl) 
31     {
32       cerr << "Can't load RunLoader from file!\n";
33       rc = 1;
34       return rc;
35     }
36
37   rl->LoadgAlice();
38   gAlice=rl->GetAliRun();
39
40   if (!gAlice)
41     {
42       cerr << "<TOFquickanal> AliRun object not found on file \n";
43       rc = 2;
44       return rc;
45     }
46
47   // Get the pointer to the TOF detector
48   AliLoader *tofl = rl->GetLoader("TOFLoader");
49   AliTOF * tof = (AliTOF*) gAlice->GetDetector("TOF");
50   if (tof == 0x0 || tofl == 0x0) {
51     cerr << "<TOFquickanal> Can not find TOF or TOFLoader\n";
52     rc = 3;
53     return rc;
54   }
55
56   //=======> Create histograms
57   //---> Time of Flight for Primary Particles (ns)
58   TH1F *htofprim = new TH1F("htofprim","Time of Flight for Primary Particles",100,0.,100.);
59   //--->Time of Flight for Secondary Particles (ns)
60   TH1F *htofsec  = new TH1F("htofsec","Time of Flight for Secondary Particles",100,0.,100.);
61   
62   //---> r (radius) coordinate of production in the ALICE frame for secondary particles that produce at 
63   //     least one TOF-hit (cm) - cylindrical coordinate system assumed, primary plus secondary-
64   TH1F *hradius = new TH1F("hradius","r (radius) coordinate at the production vertex for secondary particles with at least one TOF-Hit",50,0.,500.);
65   
66   //---> Momentum of primary particles that produce (at least) one TOF-hit when the hit
67   //     is produced (Gev/c)
68   TH1F *htofmom  = new TH1F("htofmom","Momentum of primary particles when the Hit is produced",50,0.,5.);
69   
70   //---> Momentum of primary particles that produce (at least) one TOF-hit at the production vertex
71   //     (Gev/c)
72   TH1F *hprodmom  = new TH1F("hprodmom","Momentum of primary particles (with at least one TOF hit) at the production ",50,0.,5.); 
73   
74   //---> Theta of production for primary particles that produce (at least) one TOF-hit (deg)
75   TH1F *hprodthe  = new TH1F("hprodthe","Theta of primary particles (with at least one TOF hit) at the production ",90,0.,180.);
76   
77   //---> Phi of production for primary particles that produce (at least) one TOF-hit (deg)
78   TH1F *hprodphi  = new TH1F("hprodphi","Phi of primary particles (with at least one TOF hit) at the production ",180,-180.,180.);
79   
80   //---> z Coordinate of the TOF Hit (z beam axis) - primary plus secondary - (cm)
81   TH1F *hzcoor = new TH1F("hzcoor","z Coordinate of the TOF Hit",800,-400.,400.);
82   
83   //---> Incidence Angle of the particle on the pad (or strip) (deg)  - primary plus secondary - 
84   TH1F *hincangle = new TH1F("hincangle","Incidence Angle of the particle on the strip",90,0.,180.);
85   
86   printf ("Processing event %d \n", eventNumber);
87   rl->GetEvent(eventNumber);
88   
89   // Get pointers to Alice detectors and Hits containers
90   tofl->LoadHits();
91   TTree *TH = tofl->TreeH();
92   tof->SetTreeAddress();
93   if (!TH) {
94     cout << "<TOFquickanal> No hit tree found" << endl;
95     rc = 4;
96     return rc;
97   }
98   
99   // Import the Kine Tree for the event eventNumber in the file  
100   rl->LoadHeader();
101   rl->LoadKinematics();
102   //AliStack * stack = rl->Stack();
103   
104   Int_t ntracks = TH->GetEntries();
105   cout<<" ntracks = "<<ntracks<<endl;
106   
107   AliTOFhitT0 *tofHit;
108   
109   // Start loop on tracks in the hits containers
110   for (Int_t track=0; track<ntracks;track++) {
111     
112     tof->ResetHits();
113     TH->GetEvent(track);
114     
115     for(tofHit=(AliTOFhitT0*)tof->FirstHit(track); tofHit; tofHit=(AliTOFhitT0*)tof->NextHit()) {
116       
117       Float_t toflight = tofHit->GetTof();
118       toflight        *= 1.E+09;  // conversion from s to ns
119       Double_t tofmom  = tofHit->GetMom();
120       
121       Int_t ipart = tofHit->Track();
122       TParticle *particle = gAlice->Particle(ipart);
123       if (particle->GetFirstMother() < 0) {
124         htofprim->Fill(toflight);
125         htofmom->Fill(tofmom); 
126       } else {
127         htofsec->Fill(toflight); 
128       }
129       
130       Double_t zcoor = tofHit->Z();
131       hzcoor->Fill(zcoor);
132       
133       Double_t incangle = tofHit->GetIncA();
134       hincangle->Fill(incangle);
135       
136       Double_t xcoor  = particle->Vx();
137       Double_t ycoor  = particle->Vy();
138       Double_t radius = TMath::Sqrt(xcoor*xcoor+ycoor*ycoor);
139       if (particle->GetFirstMother() >= 0) hradius->Fill(radius);
140       
141       Double_t prodmom = particle->P();        
142       if (prodmom!=0.) {
143         Double_t dummy = (particle->Pz())/prodmom;
144         Double_t prodthe = TMath::ACos(dummy);
145         prodthe *= 57.29578; // conversion from rad to deg
146         if (particle->GetFirstMother() < 0) hprodthe->Fill(prodthe);
147       } // theta at production vertex
148       
149       if (particle->GetFirstMother() < 0) {         
150         hprodmom->Fill(prodmom);
151         Double_t dummypx = particle->Px();
152         Double_t dummypy = particle->Py();
153         Double_t prodphi = TMath::ATan2(dummypy,dummypx);
154         prodphi *= 57.29578; // conversion from rad to deg
155         hprodphi->Fill(prodphi);
156       } // phi at production vertex
157     } // close loop on TOF-hits
158   } // close loop on tracks in the hits containers
159   
160   //Create  canvas, set the view range, show histograms
161   TCanvas *c1 = new TCanvas("c1","Alice TOF hits quick analysis",400,10,600,700);
162   c1->cd();
163   hprodmom->Draw();
164   
165   TCanvas *c2 = new TCanvas("c2","Alice TOF hits quick analysis",400,10,600,700);
166   c2->cd();
167   hprodthe->Draw();
168   
169   TCanvas *c3 = new TCanvas("c3","Alice TOF hits quick analysis",400,10,600,700);
170   c3->cd();
171   hprodphi->Draw();
172   
173   TCanvas *c4 = new TCanvas("c4","Alice TOF hits quick analysis",400,10,600,700);
174   c4->cd();
175   hzcoor->Draw();
176   
177   TCanvas *c5 = new TCanvas("c5","Alice TOF hits quick analysis",400,10,600,700);
178   c5->cd();
179   hradius->Draw();
180   
181   TCanvas *c6 = new TCanvas("c6","Alice TOF hits quick analysis",400,10,600,700);
182   c6->cd();
183   htofprim->Draw();
184   
185   TCanvas *c7 = new TCanvas("c7","Alice TOF hits quick analysis",400,10,600,700);
186   c7->cd();
187   htofsec->Draw();
188   
189   
190   TCanvas *c8 = new TCanvas("c8","Alice TOF hits quick analysis",400,10,600,700);
191   c8->cd();
192   htofmom->Draw();
193   
194   TCanvas *c9 = new TCanvas("c9","Alice TOF hits quick analysis",400,10,600,700);
195   c9->cd();
196   hincangle->Draw();
197   
198   //tofl->UnloadHits();
199   //rl->UnloadHeader();
200   //rl->UnloadgAlice();
201   //rl->UnloadKinematics();
202
203   return rc;
204
205 }