]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCDisplayClusters_new.C
cluster information
[u/mrichter/AliRoot.git] / TPC / AliTPCDisplayClusters_new.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
8 #endif
9 Int_t AliTPCDisplayClusters(Int_t eventn=0, Int_t noiseth=15) {
10    cerr<<"Displaying clusters...\n";
11
12    TFile *file=TFile::Open("galice.root");
13    if (!file->IsOpen()) {cerr<<"Can't open galice.root !\n"; return 1;}
14
15    TFile *cf=TFile::Open("AliTPCclusters.root");
16    if (!cf->IsOpen()){cerr<<"Can't open AliTPCclusters.root !\n"; return 3;}
17
18    AliTPCParam *dig=(AliTPCParam *)cf->Get("75x40_100x60");
19    if (!dig) {cerr<<"TPC parameters have not been found !\n"; return 2;}
20
21
22    // some "constants"
23    Int_t markerColorSignal = 5;
24    Int_t markerColorBgr = 2;
25    Int_t MASK = 10000000;
26
27
28    TCanvas *c1=new TCanvas("cdisplay", "Cluster display",0,0,700,730);
29    TView *v=new TView(1);
30    v->SetRange(-430,-560,-430,430,560,1710);
31    c1->Clear();
32    c1->SetFillColor(1);
33    c1->SetTheta(90.);
34    c1->SetPhi(0.);
35
36    AliTPCClustersArray *ca=new AliTPCClustersArray;
37    ca->Setup(dig);
38    ca->SetClusterType("AliTPCcluster");
39    char  cname[100];
40    sprintf(cname,"TreeC_TPC_%d",eventn);
41
42    ca->ConnectTree(cname);
43    Int_t nrows=Int_t(ca->GetTree()->GetEntries());
44    for (Int_t n=0; n<nrows; n++) {
45        AliSegmentID *s=ca->LoadEntry(n);
46        Int_t sec,row;
47        dig->AdjustSectorRow(s->GetID(),sec,row);
48        AliTPCClustersRow &clrow = *ca->GetRow(sec,row);
49        Int_t ncl=clrow.GetArray()->GetEntriesFast();
50        TPolyMarker3D *pm=new TPolyMarker3D(ncl);
51        TPolyMarker3D *pmSignal=new TPolyMarker3D(ncl); // polymarker for signal
52        Int_t imarBgr=0;
53        Int_t imarSignal=0;
54        while (ncl--) {
55            AliTPCcluster *cl=(AliTPCcluster*)clrow[ncl];
56            Double_t x=dig->GetPadRowRadii(sec,row), y=cl->GetY(), z=cl->GetZ();
57            if (cl->GetQ()<noiseth) continue;
58            Float_t cs, sn, tmp;
59            dig->AdjustCosSin(sec,cs,sn);
60            tmp = x*cs-y*sn; y= x*sn+y*cs; x=tmp; 
61            if (TMath::Abs(z)>50) continue;
62            Int_t trackId = cl->GetLabel(0);
63            if (trackId<0) continue;  //bug in tracks ID
64            if (trackId < MASK-1) {
65              pmSignal->SetPoint(imarSignal,x,y,z);
66              imarSignal++;
67            } else {
68              pm->SetPoint(imarBgr,x,y,z);
69              imarBgr++;
70            }          
71        }
72        ca->ClearRow(sec,row);
73        
74       // change color for signal
75        pm->SetMarkerSize(1); 
76        pm->SetMarkerColor(markerColorBgr);
77        pm->SetMarkerStyle(1);
78        pm->Draw();
79
80        pmSignal->SetMarkerSize(1); 
81        pmSignal->SetMarkerColor(markerColorSignal);
82        pmSignal->SetMarkerStyle(1);
83        pmSignal->Draw();      
84    }
85    delete ca;
86    cf->Close();
87
88    TGeometry *geom=(TGeometry*)file->Get("AliceGeom");
89    //   TList *list = geom->GetListOfNodes();
90    TNode * main = (TNode*)((geom->GetListOfNodes())->First());
91    TIter next(main->GetListOfNodes());
92    TNode  *module=0;
93    while((module = (TNode*)next())) {
94      char ch[100];
95      sprintf(ch,"%s\n",module->GetTitle());
96      //printf("%s\n",module->GetTitle());
97      if (ch[0]=='T'&&ch[1]=='P' && ch[2]=='C')  //if TPC draw
98        module->SetVisibility(3);
99      else
100        module->SetVisibility(-1);
101    }
102      
103    
104    geom->Draw("same");
105    //v->Draw();
106    c1->Modified(); c1->Update(); 
107
108    file->Close();
109    return 0;
110 }
111