]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCDisplayClusters.C
Removing old classes from LinkDef
[u/mrichter/AliRoot.git] / TPC / AliTPCDisplayClusters.C
1 /****************************************************************************
2  *           Origin: I.Belikov, CERN, Jouri.Belikov@cern.ch                 *
3  ****************************************************************************/
4
5 #ifndef __CINT__
6 #include <iostream.h>
7 #include <TFile.h>
8 #include <TTree.h>
9 #include <TCanvas.h>
10 #include <TView.h>
11 #include <TPolyMarker3D.h>
12 #include <TNode.h>
13 #include <TGeometry.h>
14
15 #include "AliTPCParam.h"
16 #include "AliTPCClustersArray.h"
17 #include "AliTPCClustersRow.h"
18 #include "AliTPCcluster.h"
19 #endif
20
21 Int_t AliTPCDisplayClusters(Int_t eventn=0, Int_t noiseth=15) {
22    cerr<<"Displaying clusters...\n";
23
24    TFile *file=TFile::Open("galice.root");
25    if (!file->IsOpen()) {cerr<<"Can't open galice.root !\n"; return 1;}
26
27    TFile *cf=TFile::Open("AliTPCclusters.root");
28    if (!cf->IsOpen()){cerr<<"Can't open AliTPCclusters.root !\n"; return 3;}
29
30    AliTPCParam *dig=(AliTPCParam *)cf->Get("75x40_100x60_150x60");
31    if (!dig) {cerr<<"TPC parameters have not been found !\n"; return 2;}
32
33    TCanvas *c1=new TCanvas("cdisplay", "Cluster display",0,0,700,730);
34    TView *v=new TView(1);
35    v->SetRange(-430,-560,-430,430,560,1710);
36    c1->Clear();
37    c1->SetFillColor(1);
38    c1->SetTheta(90.);
39    c1->SetPhi(0.);
40
41    AliTPCClustersArray *ca=new AliTPCClustersArray;
42    ca->Setup(dig);
43    ca->SetClusterType("AliTPCcluster");
44    char  cname[100];
45    sprintf(cname,"TreeC_TPC_%d",eventn);
46
47    ca->ConnectTree(cname);
48    Int_t nrows=Int_t(ca->GetTree()->GetEntries());
49    for (Int_t n=0; n<nrows; n++) {
50        AliSegmentID *s=ca->LoadEntry(n);
51        Int_t sec,row;
52        dig->AdjustSectorRow(s->GetID(),sec,row);
53        AliTPCClustersRow &clrow = *ca->GetRow(sec,row);
54        Int_t ncl=clrow.GetArray()->GetEntriesFast();
55        TPolyMarker3D *pm=new TPolyMarker3D(ncl);
56        while (ncl--) {
57            AliTPCcluster *cl=(AliTPCcluster*)clrow[ncl];
58            Double_t x=dig->GetPadRowRadii(sec,row), y=cl->GetY(), z=cl->GetZ();
59            if (cl->GetQ()<noiseth) continue;
60            Float_t cs, sn, tmp;
61            dig->AdjustCosSin(sec,cs,sn);
62            tmp = x*cs-y*sn; y= x*sn+y*cs; x=tmp;
63            pm->SetPoint(ncl,x,y,z);
64        }
65        ca->ClearRow(sec,row);
66        pm->SetMarkerSize(1); pm->SetMarkerColor(2); pm->SetMarkerStyle(1);
67        pm->Draw();
68    }
69    delete ca;
70    cf->Close();
71
72    TGeometry *geom=(TGeometry*)file->Get("AliceGeom");
73    TNode * main = (TNode*)((geom->GetListOfNodes())->First());
74    TIter next(main->GetListOfNodes());
75    TNode  *module=0;
76    while((module = (TNode*)next())) {
77      char ch[100];
78      sprintf(ch,"%s\n",module->GetTitle());
79      //printf("%s\n",module->GetTitle());
80      if (ch[0]=='T'&&ch[1]=='P' && ch[2]=='C')  //if TPC draw
81        module->SetVisibility(3);
82      else
83        module->SetVisibility(-1);
84    }
85      
86    
87    geom->Draw("same");
88    //v->Draw();
89    c1->Modified(); c1->Update(); 
90
91    file->Close();
92    return 0;
93 }
94