]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCDisplayDigits.C
Generation of dictionaries and rootmaps with Root6
[u/mrichter/AliRoot.git] / TPC / AliTPCDisplayDigits.C
1 /// \file AliTPCDisplayDigits.C
2 /// \author I.Belikov, CERN, Jouri.Belikov@cern.ch
3
4 #ifndef __CINT__
5 #include <Riostream.h>
6 #include <TFile.h>
7 #include <TTree.h>
8 #include <TCanvas.h>
9 #include <TStyle.h>
10 #include <TH2.h>
11
12 #include "AliTPCParam.h"
13 #include "AliSimDigits.h"
14 #endif
15
16 Int_t AliTPCDisplayDigits(Int_t eventn, int sec, int row, int lab=-1,
17                  int max_t_chan=500, float min_t=0., float max_t=500.,
18                  int max_p_chan=150, float min_p=0., float max_p=150.)
19 {
20    cerr<<"Displaying digits...\n";
21
22 // Connect the Root Galice file containing Geometry, Kine and Hits
23    TFile *f = TFile::Open("rfio:galice.root");
24    if (!f->IsOpen()) {cerr<<"Can't open galice.root !\n"; return 1;}
25
26    
27    
28    AliTPCParam *par=(AliTPCParam *)f->Get("75x40_100x60_150x60");
29    if (!par) { cerr<<"TPC parameters have not been found !\n"; return 2; }
30
31    Char_t s[80];
32    sprintf(s,"Sector %d   Row %d\n",sec,row);
33    TH2F *h = new TH2F("h",s,max_t_chan,min_t,max_t,max_p_chan,min_p,max_p);
34    
35    Char_t  name[100];
36    sprintf(name,"TreeD_75x40_100x60_150x60_%d",eventn);
37
38    TTree *t=(TTree*)f->Get(name);
39    if (!t) { cerr<<"TPC digits have not been found !\n"; return 3; }
40    AliSimDigits dummy, *digit=&dummy;
41    t->GetBranch("Segment")->SetAddress(&digit);
42    Int_t sbr=(Int_t)t->GetEntries();
43    for (Int_t i=0; i<sbr; i++) {
44        if (!t->GetEvent(i)) continue;
45        Int_t s,r;
46        par->AdjustSectorRow(digit->GetID(),s,r);
47        if (s==sec && r==row) goto ok;
48    }
49    return 4;
50    
51 ok:
52
53    Int_t imax=0, jmax=0, qmax=0;
54    digit->First();
55    do {
56      //Short_t dig=digit->CurrentDigit();
57       Int_t i=digit->CurrentRow(), j=digit->CurrentColumn();
58       Short_t dig=digit->GetDigit(i,j);
59       if (lab >= 0) {
60          Int_t lab0=digit->GetTrackID(i,j,0);
61          Int_t lab1=digit->GetTrackID(i,j,1);
62          Int_t lab2=digit->GetTrackID(i,j,2);
63          if (lab0!=lab) if (lab1!=lab) if (lab2!=lab) continue;
64          if (dig>qmax) {imax=i; jmax=j; qmax=dig;}
65          cerr<<lab0<<' '<<lab1<<' '<<lab2<<endl;
66       }
67       h->Fill(i,j,dig);
68    } while (digit->Next());
69    if (qmax>0){cerr<<"Peak (time,pad,q) : "<<imax<<' '<<jmax<<' '<<qmax<<endl;}
70
71    h->SetMaximum(100);
72    gStyle->SetOptStat(0);
73    TCanvas *c1=new TCanvas("ddisplay","TPC digits display",0,0,1110,680);
74    TPad *p1=new TPad("p1","",0,0,1,0.5);
75    p1->Draw();
76    TPad *p2=new TPad("p2","",0,0.5,1,1);
77    p2->Draw();
78    p2->cd();
79    h->DrawCopy("lego");
80    p1->cd();
81    h->DrawCopy("colz");
82
83    c1->Modified(); c1->Update(); 
84
85    f->Close();
86    return 0;
87 }
88