]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCDisplayDigits3Dnew.C
GetClusterFast function implemented (No getter before) (Marian)
[u/mrichter/AliRoot.git] / TPC / AliTPCDisplayDigits3Dnew.C
1 #ifndef __CINT__
2   #include "alles.h"
3   #include "AliTPCtracker.h"
4   #include "TView.h"
5   #include "TPolyMarker3D.h"
6   #include "AliSimDigits.h"
7   #include "AliRunLoader.h"
8   #include "AliLoader.h"
9   #include "AliTPCParamSR.h"
10 #endif
11
12 //  
13 //  display 3D digits
14 //  input parameter is event number  - threshol to the noise
15 //  if sdigits=kTRUE it will display Sdigits - otherwise it display digits
16 //  signal event is displayed with yellow color
17
18 Int_t AliTPCDisplayDigits3Dnew(Int_t eventn=0, Int_t noiseth=15, Bool_t sdigits=kFALSE) {
19    cerr<<"Displaying digits...\n";
20
21 //   TFile *file=TFile::Open("galice.root");
22 //   if (!file->IsOpen()) {cerr<<"Can't open galice.root !\n"; return 1;}
23
24 //   TFile *cf=TFile::Open("galice.root");
25 //    if (!cf->IsOpen()){cerr<<"Can't open AliTPCclusters.root !\n"; return 3;}
26
27    AliRunLoader* rl = AliRunLoader::Open();
28    rl->GetEvent(eventn);
29
30    AliLoader* tpcl = (AliLoader*)rl->GetLoader("TPCLoader");
31    if (tpcl == 0x0)
32     {
33       cerr<<"Can not get TPC Loader"<<endl;
34       return 1;
35     }
36
37    rl->CdGAFile();
38    AliTPCParam *param=(AliTPCParam *)gDirectory->Get("75x40_100x60");
39    if(param){
40      cerr<<"2 pad-length geom hits with 3 pad-lengths geom parameters\n";
41      delete param;
42      param = new AliTPCParamSR();
43    }
44    else
45    {
46      param=(AliTPCParamSR *)gDirectory->Get("75x40_100x60_150x60");
47    }
48
49    if (!param) {cerr<<"TPC parameters have not been found !\n"; return 2;}
50
51    TCanvas *c1=new TCanvas("ddisplay", "Digits display",0,0,700,730);
52    TView *v=new TView(1);
53    v->SetRange(-430,-560,-430,430,560,1710);
54    c1->Clear();
55    c1->SetFillColor(1);
56    c1->SetTheta(90.);
57    c1->SetPhi(0.);
58
59    AliTPCDigitsArray *digarr=new AliTPCDigitsArray;
60    digarr->Setup(param);
61    
62   
63   TTree* tree;
64    
65    char  cname[100];
66    if (!sdigits)
67      {
68        tpcl->LoadDigits();
69        tree = tpcl->TreeD();
70      }
71    else
72      {
73        tpcl->LoadSDigits();
74        tree = tpcl->TreeS();
75      }
76
77 // some "constants"
78    Int_t markerColorSignal = 5;
79    Int_t markerColorBgr = 2;
80    Int_t MASK = 10000000;
81
82    Int_t imarBgr=0;
83    Int_t imarSignal=0;
84    Int_t ncl=0;
85    
86    digarr->ConnectTree(tree);
87    
88    Int_t nrows=Int_t(digarr->GetTree()->GetEntries());
89    Int_t all0=0;
90    for (Int_t n=0; n<nrows; n++) {
91        AliSimDigits *s=(AliSimDigits*)digarr->LoadEntry(n);
92        Int_t sec,row;
93        param->AdjustSectorRow(s->GetID(),sec,row);
94        Int_t npads, sign;
95        {
96        const Int_t kNIS=param->GetNInnerSector(), kNOS=param->GetNOuterSector();
97        if (sec < kNIS) {
98          npads = param->GetNPadsLow(row);
99          sign = (sec < kNIS/2) ? 1 : -1;
100        } else {
101          npads = param->GetNPadsUp(row);
102          sign = ((sec-kNIS) < kNOS/2) ? 1 : -1;
103        }
104        }
105
106        AliSimDigits *digrow = (AliSimDigits*)digarr->GetRow(sec,row);
107        ncl=0;
108        if (digrow->First()){
109        while(digrow->Next()) ncl++;
110        ncl++;
111        }
112        TPolyMarker3D *pm=new TPolyMarker3D(ncl);
113        TPolyMarker3D *pmSignal=new TPolyMarker3D(ncl); // polymarker for signal
114        imarBgr=0;
115        imarSignal=0;
116        if (digrow->First()) do {
117        Short_t dig=digrow->CurrentDigit();
118        Double_t y = (digrow->CurrentColumn()- 0.5 - 0.5*npads)*param->GetPadPitchWidth(sec);
119        Double_t z = sign*(param->GetZLength()-param->GetZWidth()*digrow->CurrentRow());       
120        Double_t x=param->GetPadRowRadii(sec,row);
121 //       cout<<dig<<endl;
122        if (dig<noiseth) continue;
123 //       cout<<"\nAbove noise Threshold";
124        Int_t trackId = digrow->GetTrackID(digrow->CurrentRow(),digrow->CurrentColumn(),0);      
125        Float_t cs, sn, tmp;
126        param->AdjustCosSin(sec,cs,sn);
127        tmp = x*cs-y*sn; y= x*sn+y*cs; x=tmp;
128        if (trackId<0) continue;  //if signal from croostalk - we don't have track ID information
129 //       cout<<"Track ID > 0";
130              if (trackId < MASK-1) {
131          pmSignal->SetPoint(imarSignal,x,y,z);
132          imarSignal++;
133        } else {
134          pm->SetPoint(imarBgr,x,y,z);
135          imarBgr++;
136        }
137        } while (digrow->Next());
138
139 // change color for signal
140        digarr->ClearRow(sec,row);
141        pm->SetMarkerSize(1); 
142        pm->SetMarkerColor(markerColorBgr);
143        pm->SetMarkerStyle(1);
144        pm->Draw();
145
146        pmSignal->SetMarkerSize(1); 
147        pmSignal->SetMarkerColor(markerColorSignal);
148        pmSignal->SetMarkerStyle(1);
149        pmSignal->Draw();
150        all0+=imarSignal;
151 //       cout<<"imarSignal ="<<imarSignal<<"   imarBgr ="<<imarBgr<<"   ncl ="<<ncl<<endl;
152    }
153    printf("%d\n",all0);
154    
155    
156
157    delete digarr;
158
159    //draw TPC skeleton
160    rl->CdGAFile();
161    TGeometry *geom=(TGeometry*)gDirectory->Get("AliceGeom");
162    TNode * main = (TNode*)((geom->GetListOfNodes())->First());
163    TIter next(main->GetListOfNodes());
164    TNode  *module=0;
165    while((module = (TNode*)next())) {
166      char ch[100];
167      sprintf(ch,"%s\n",module->GetTitle());
168      //printf("%s\n",module->GetTitle());
169      if (ch[0]=='T'&&ch[1]=='P' && ch[2]=='C')  //if TPC draw
170        module->SetVisibility(3);
171      else
172        module->SetVisibility(-1);
173    }
174      
175    
176    geom->Draw("same");
177    //v->Draw();
178    c1->Modified(); c1->Update(); 
179    
180    delete rl;
181    return 0;
182 }
183