]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCDisplayClustersMI.C
adding new library libAliHLTGlobal for global HLT components; adding global track...
[u/mrichter/AliRoot.git] / TPC / AliTPCDisplayClustersMI.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include "alles.h"
3 #include "AliTPCtracker.h"
4 #include "TView3D.h"
5 #include "TPolyMarker3D.h"
6 #include "AliSimDigits.h"
7 #include "AliTPCParam.h"
8 #include "AliRunLoader.h"
9 #include "AliTPCclusterMI.h"
10 #endif
11
12 /*
13   Author:   marian.ivanov@cern.ch
14
15   How to use ?
16   TGeoManager::Import("geometry.root");
17   .L AliTPCDisplayClustersMI.C+
18   AliTPCDisplayClusters disp;
19   disp.Init(0,0);   //specify event number, and threshold for the noise
20   disp.DisplayClusters();  
21
22 */
23
24
25
26
27 class AliTPCDisplayClusters{
28 public:
29   AliTPCDisplayClusters();
30   void SetIO(Int_t event);
31   void LoadClusters(Int_t noiseth);  
32   void DisplayClusters(Int_t first=0, Int_t last=-1);
33   void Init(Int_t event, Int_t noiseth){SetIO(event); LoadClusters(noiseth);}
34   TObjArray   * fArray;
35   AliTPCParam * fParam;
36   TTree       * fTree;
37   TGeometry   * fGeom;
38 };
39
40
41 //----------------------------------------------------------------------
42 AliTPCDisplayClusters::AliTPCDisplayClusters()
43 {
44   fArray = 0;
45   fParam = 0;
46   fTree  = 0;
47   fGeom  = 0;
48 }
49 //----------------------------------------------------------------------
50 void AliTPCDisplayClusters::SetIO(Int_t event)
51 {
52   AliRunLoader* rl = AliRunLoader::Open();
53   rl->GetEvent(event);  
54   AliLoader* tpcl = (AliLoader*)rl->GetLoader("TPCLoader");
55   if (tpcl == 0x0)
56     {
57       cerr<<"Can not get TPC Loader"<<endl;
58       return;
59     }  
60   rl->CdGAFile();
61   fParam=(AliTPCParam *)gDirectory->Get("75x40_100x60_150x60");
62   fGeom=(TGeometry*)gDirectory->Get("AliceGeom");
63   if(!fParam){
64     fParam = new AliTPCParamSR();
65     fParam->Update();
66   }
67   
68   if (!fParam) {cerr<<"TPC parameters have not been found !\n"; return ;}  
69   tpcl->LoadRecPoints();
70   fTree = tpcl->TreeR();
71    
72 }
73 //----------------------------------------------------------------------
74 void AliTPCDisplayClusters::LoadClusters(Int_t noiseth)
75 {
76   //
77   // load all clusters to memory
78   if (fArray) {
79     fArray->Delete();
80     delete fArray;
81   }
82   fArray = new TObjArray(fParam->GetNSegmentsTotal());
83   //
84   AliTPCClustersRow * pclrow = new AliTPCClustersRow;
85   pclrow->SetClass("AliTPCclusterMI");
86   pclrow->SetArray(0);
87   TBranch *br = fTree->GetBranch("Segment");
88   br->SetAddress(&pclrow);
89   //
90   //
91   Int_t nrows=Int_t(fTree->GetEntries());
92   for (Int_t n=0; n<nrows; n++) {
93     //
94     pclrow = new AliTPCClustersRow;
95     pclrow->SetClass("AliTPCclusterMI");
96     pclrow->SetArray(0);
97     br->SetAddress(&pclrow);
98     //
99     br->GetEntry(n);
100     AliTPCClustersRow &clrow = *pclrow;
101     Int_t ncl=clrow.GetArray()->GetEntriesFast();
102     TObjArray * arrrow = new TObjArray(0);
103     fArray->AddAt(arrrow,pclrow->GetID());
104
105     // printf("segment\t%d\trow\t%d\tclusters\t%d",n,pclrow->GetID(),ncl);
106     while (ncl--) {
107       AliTPCclusterMI *cl=(AliTPCclusterMI*)clrow[ncl];
108       if (cl->GetQ()>noiseth){
109         arrrow->AddLast(cl);    
110       }  
111     }   
112     //    printf("over\t%d\n",arrrow->GetEntries());
113
114   }
115 }
116
117
118 void AliTPCDisplayClusters::DisplayClusters(Int_t first, Int_t last)
119 {
120   Int_t nrows = fParam->GetNSegmentsTotal();
121   // some "constants"
122   Int_t markerColorSignal = 5;
123   Int_t markerColorBgr = 2;
124   Int_t MASK = 10000000;
125
126   TCanvas *c1=new TCanvas("cdisplay", "Cluster display",0,0,700,730);
127   TView3D *v=new TView3D();
128   v->SetRange(-330,-360,-330,360,360,1710);
129   c1->Clear();
130   c1->SetFillColor(1);
131   c1->SetTheta(90.);
132   c1->SetPhi(0.);
133   
134   for (Int_t irow=0; irow<nrows; irow++) {
135
136     TObjArray * arr = (TObjArray*)fArray->At(irow);
137     if (!arr) continue;
138     
139     Int_t sec,row;
140     fParam->AdjustSectorRow(irow,sec,row);    
141     Int_t ncl=arr->GetEntriesFast();
142
143     TPolyMarker3D *pm=new TPolyMarker3D(ncl);
144     TPolyMarker3D *pmSignal=new TPolyMarker3D(ncl); // polymarker for signal
145     Int_t imarBgr=0;
146     Int_t imarSignal=0;
147     
148     while (ncl--) {
149       AliTPCclusterMI *cl=(AliTPCclusterMI*)arr->At(ncl);
150       if (cl){
151       //
152         Double_t x=fParam->GetPadRowRadii(sec,row), y=cl->GetY(), z=cl->GetZ();
153         Float_t cs, sn, tmp;
154         fParam->AdjustCosSin(sec,cs,sn);
155         tmp = x*cs-y*sn; y= x*sn+y*cs; x=tmp; 
156         Float_t clxyz[3];
157         cl->GetGlobalXYZ(clxyz);
158
159         Int_t trackId = cl->GetLabel(0);
160         if ( (last>0) &&trackId>last) continue;
161         if (trackId<first) continue;
162         if (trackId < MASK-1) {
163           pmSignal->SetPoint(imarSignal,clxyz[0],clxyz[1],clxyz[2]);
164           imarSignal++;
165         } else {
166           pm->SetPoint(imarBgr,clxyz[0],clxyz[1],clxyz[2]);
167           imarBgr++;
168         }          
169       }
170     }
171         
172     // change color for signal
173     pm->SetMarkerSize(1); 
174     pm->SetMarkerColor(markerColorBgr);
175     pm->SetMarkerStyle(1);
176     pm->Draw();
177     
178     pmSignal->SetMarkerSize(1); 
179     pmSignal->SetMarkerColor(markerColorSignal);
180     pmSignal->SetMarkerStyle(1);
181     pmSignal->Draw();      
182   }
183   
184   
185  //  TNode * main = (TNode*)((fGeom->GetListOfNodes())->First());
186 //   TIter next(main->GetListOfNodes());
187 //   TNode  *module=0;
188 //   while((module = (TNode*)next())) {
189 //     char ch[100];
190 //     sprintf(ch,"%s\n",module->GetTitle());
191 //     //printf("%s\n",module->GetTitle());
192 //     if (ch[0]=='T'&&ch[1]=='P' && ch[2]=='C')  //if TPC draw
193 //       module->SetVisibility(3);
194 //     else
195 //       module->SetVisibility(-1);
196 //   }
197   
198   
199 //   fGeom->Draw("same");
200   c1->Modified(); c1->Update(); 
201   
202   
203
204 }