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