]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/analITS.C
Delete only existent objects.
[u/mrichter/AliRoot.git] / ITS / analITS.C
1 void analITS (Int_t evNumber=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.
7 //   
8 //     Root > .L analITS.C   //this loads the macro in memory
9 //     Root > analITS();     //by default process first event   
10 //     Root > analITS(2);    //process third event
11 //Begin_Html
12 /*
13 <img src="figures/analITS_ref.gif">
14 </pre>
15 <br clear=left>
16 <font size=+2 color=red>
17 <p>A reference plot produced by analITS.C.
18 </font>
19 <pre>
20  */
21 //End_Html
22 /////////////////////////////////////////////////////////////////////////
23
24
25 // Dynamically link some shared libs
26    if (gClassTable->GetID("AliRun") < 0) {
27       gROOT->LoadMacro("loadlibs.C");
28       loadlibs();
29    } // end if gClassTable...
30       
31 // Connect the Root Galice file containing Geometry, Kine and Hits
32    TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject("galice.root");
33    if (!file) file = new TFile("galice.root");
34
35 // Get AliRun object from file or create it if not on file
36    if (!gAlice) {
37       gAlice = (AliRun*)file->Get("gAlice");
38       if (gAlice) printf("AliRun object found on file\n");
39       if (!gAlice) gAlice = new AliRun("gAlice","Alice test program");
40    }
41       
42 // Import the Kine and Hits Trees for the event evNumber in the file
43    Int_t nparticles = gAlice->GetEvent(evNumber);
44    if (nparticles <= 0) return;
45    Float_t x,y,z,mass,e,r;
46    Int_t nbytes = 0;
47    Int_t j,hit,ipart;
48    Int_t nhits;
49    Int_t sector,plane;
50    TParticle *particle;
51    AliTPChit  *tpcHit;
52    AliITShit  *itsHit;
53
54 // Get pointers to Alice detectors and Hits containers
55    AliDetector *TPC  = gAlice->GetDetector("TPC");
56    AliDetector *ITS  = gAlice->GetDetector("ITS");
57    TClonesArray *Particles = gAlice->Particles();
58    if (TPC) TClonesArray *TPChits   = TPC->Hits();
59    if (ITS) TClonesArray *ITShits   = ITS->Hits();
60
61    TTree *TH = gAlice->TreeH();
62    Int_t ntracks    = TH->GetEntries();
63
64    // Create histograms
65    TH1F *hITSZ    = new TH1F("hITSZ","Z of hits in ITS",100,-50.,50.);
66    TH1F *hITSR    = new TH1F("hITSR","R of hits in ITS",100,0.,50.);
67    TH1F *hITSlayer= new TH1F("hITSlayer","JLAY of hits in ITS",6, 0., 6.);
68    TH1F *hITSTr   = new TH1F("hITSTr","Track number for hits in ITS",100,0.,50000.);
69
70 // Start loop on tracks in the hits containers
71    for (Int_t track=0; track<ntracks;track++) {
72      gAlice->ResetHits();
73      nbytes += TH->GetEvent(track);
74      Int_t i=0;
75      if (ITS) {
76        nhits = ITShits->GetEntriesFast();
77        for (hit=0;hit<nhits;hit++) {
78          itsHit   = (AliITShit*)ITShits->UncheckedAt(hit);
79          // With this new version, to be able to do proper detector
80          // simulations, the requirment that a "hit" leave some
81          // energy deposited has been removed.
82          if(itsHit->GetIonization()<=0.0) continue;
83          ipart    = itsHit->fTrack;
84          particle = (TParticle*)Particles->UncheckedAt(ipart);
85          z = itsHit->GetZG();
86          hITSZ->Fill(z);
87          r = sqrt(itsHit->GetXG()*itsHit->GetXG() + 
88                   itsHit->GetYG()*itsHit->GetYG());
89          hITSR->Fill(r);
90          hITSlayer->Fill((Float_t)itsHit->GetLayer());
91          hITSTr->Fill(itsHit->fTrack);
92          i++;
93        }
94      }        
95    }
96
97 //Create a canvas, set the view range, show histograms
98    TCanvas *c1 = new TCanvas("c1","Alice TPC and ITS hits",400,10,600,700);
99    c1->Divide(2,2);
100    c1->cd(1);
101    gPad->SetFillColor(33);
102    hITSZ->SetFillColor(33);
103    hITSZ->Draw();
104    c1->cd(2);
105    gPad->SetFillColor(33);
106    hITSR->SetFillColor(46);
107    hITSR->Draw();
108    c1->cd(3);
109    gPad->SetFillColor(33);
110    hITSlayer->SetFillColor(46);
111    hITSlayer->Draw();
112    c1->cd(4);
113    gPad->SetFillColor(33);
114    hITSTr->SetFillColor(46);
115    hITSTr->Draw();
116    c1->Print("analITS.gif");
117 }