]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDdisplayDigits3D.C
Added preprocessor conditionals to support ROOT > 5.11.2.
[u/mrichter/AliRoot.git] / TRD / AliTRDdisplayDigits3D.C
1 //_____________________________________________________________________________
2 Int_t AliTRDdisplayDigits3D(Int_t event = 0, Int_t thresh = 2
3                           , Bool_t sdigits = kFALSE) 
4 {
5   //  
6   //  TRD digits display
7   //
8   //  Input parameter:
9   //    <event>   : Event number 
10   //    <thresh>  : Threshold to suppress the noise
11   //    <sdigits> : If kTRUE it will display summable digits, normal digits otherwise.
12   //                The signal event is displayed in yellow.
13   //
14
15   Char_t *inputFile = "galice.root";
16
17   const Int_t kNdict = 3;
18
19   // Define the objects
20   AliTRDv1       *trd;
21   AliTRDgeometry *geo;
22
23   Int_t           track;
24   Int_t           idict;
25
26   TString evfoldname = AliConfig::GetDefaultEventFolderName();
27   fRunLoader = AliRunLoader::GetRunLoader(evfoldname);
28   if (!fRunLoader) {
29     fRunLoader = AliRunLoader::Open(inputFile
30                                    ,AliConfig::GetDefaultEventFolderName()
31                                    ,"UPDATE");
32   }
33   if (!fRunLoader) {
34     Printf("Can not open session for file %s.",inputFile);
35     return kFALSE;
36   }
37    
38   if (!fRunLoader->GetAliRun()) {
39     fRunLoader->LoadgAlice();
40   }
41   gAlice = fRunLoader->GetAliRun();  
42   if (!gAlice) {
43     printf("Could not find AliRun object.\n");
44     return kFALSE;
45   }
46
47   fRunLoader->GetEvent(event);
48   
49   AliLoader *loader = fRunLoader->GetLoader("TRDLoader");
50   if (!loader) {
51     printf("Can not get TRD loader from Run Loader");
52   }
53   loader->LoadDigits();
54   
55   // Get the pointer to the detector object
56   trd = (AliTRDv1*) gAlice->GetDetector("TRD");
57
58   // Get the pointer to the geometry object
59   if (trd) {
60     geo = trd->GetGeometry();
61   }
62   else {
63     printf("Cannot find the geometry\n");
64     return 1;
65   }
66
67   TCanvas *c1 = new TCanvas("digits","TRD digits display",0,0,700,730);
68   TView   *v  = new TView(1);
69   v->SetRange(-430,-560,-430,430,560,1710);
70   c1->Clear();
71   c1->SetFillColor(1);
72   c1->SetTheta(90.0);
73   c1->SetPhi(0.0);
74
75   Int_t markerColorSignal = 2;
76   Int_t markerColorBgnd   = 7;
77   Int_t markerColorMerged = 5;
78   Int_t mask              = 10000000;
79
80   // Create the digits manager
81   AliTRDdigitsManager *digitsManager = new AliTRDdigitsManager();
82   digitsManager->SetDebug(1);
83   digitsManager->SetSDigits(sdigits);
84
85   // Read the digits from the file
86   if (sdigits) {
87     digitsManager->ReadDigits(loader->TreeS());
88   }
89   else {
90     if (!loader->TreeD()) {
91       printf("mist\n");
92       return kFALSE;
93     }
94     digitsManager->ReadDigits(loader->TreeD());
95   }
96
97   Int_t totalsignal = 0;
98   Int_t totalbgnd   = 0;
99   Int_t totalmerged = 0;
100
101   AliTRDparameter *par = new AliTRDparameter("TRDparameter","TRD parameter class");
102  
103   // Loop through all detectors
104   for (Int_t idet = 0; idet < geo->Ndet(); idet++) {
105
106     printf("<AliTRDdisplayDigits3D> Loading detector %d\n",idet);
107     AliTRDdataArrayI *digits  = digitsManager->GetDigits(idet);
108     digits->Expand();
109     AliTRDdataArrayI *tracks[kNdict];
110     for (Int_t idict = 0; idict < kNdict; idict++) {
111       tracks[idict] = digitsManager->GetDictionary(idet,idict);
112       tracks[idict]->Expand();
113     }
114
115     Int_t isec    = geo->GetSector(idet);
116     Int_t icha    = geo->GetChamber(idet);
117     Int_t ipla    = geo->GetPlane(idet);
118     Int_t  rowMax = par->GetRowMax(ipla,icha,isec);
119     Int_t  colMax = par->GetColMax(ipla);
120     Int_t timeMax = par->GetTimeMax();
121
122     Int_t ndigits = digits->GetOverThreshold(thresh);
123
124     if (ndigits > 0) {
125
126       TPolyMarker3D *pmSignal = new TPolyMarker3D(ndigits);
127       TPolyMarker3D *pmBgnd   = new TPolyMarker3D(ndigits);
128       TPolyMarker3D *pmMerged = new TPolyMarker3D(ndigits);
129  
130       Int_t ibgnd   = 0;
131       Int_t isignal = 0;
132       Int_t imerged = 0;
133
134       for (Int_t time = 0; time < timeMax; time++) {
135         for (Int_t  col = 0;  col <  colMax;  col++) {
136           for (Int_t  row = 0;  row <  rowMax;  row++) {
137
138             Int_t type = 1;
139
140             Int_t amp = digits->GetDataUnchecked(row,col,time);
141             for (idict = 0; idict < kNdict; idict++) {
142               Int_t trk = tracks[idict]->GetDataUnchecked(row,col,time) - 1;
143               if ((idict == 0) && (trk >= mask)) {
144                 type = 2;
145               }
146               if ((type  == 1) && (trk >= mask)) {
147                 type = 3;
148               }              
149             }
150
151             if (amp > thresh) {
152           
153               Double_t glb[3];
154               Double_t loc[3];
155
156               loc[0] = row;
157               loc[1] = col;
158               loc[2] = time;
159               geo->Local2Global(idet,loc,glb,par);
160               Double_t x = glb[0];
161               Double_t y = glb[1];
162               Double_t z = glb[2];
163
164               if      (type == 1) {
165                 pmSignal->SetPoint(isignal,x,y,z);
166                 isignal++;
167                 totalsignal++;
168               }
169               else if (type == 2) {
170                 pmBgnd->SetPoint(ibgnd,x,y,z);
171                 ibgnd++;
172                 totalbgnd++;
173               }
174               else if (type == 3) {
175                 pmMerged->SetPoint(imerged,x,y,z);
176                 imerged++;
177                 totalmerged++;
178               }
179
180             }
181
182           }
183         }
184       }
185
186       digits->Compress(1,0);
187       for (idict = 0; idict < kNdict; idict++) {
188         tracks[idict]->Compress(1,0);
189       }
190
191       pmMerged->SetMarkerSize(1); 
192       pmMerged->SetMarkerColor(markerColorMerged);
193       pmMerged->SetMarkerStyle(1);
194       pmMerged->Draw();
195
196       pmBgnd->SetMarkerSize(1); 
197       pmBgnd->SetMarkerColor(markerColorBgnd);
198       pmBgnd->SetMarkerStyle(1);
199       pmBgnd->Draw();
200
201       pmSignal->SetMarkerSize(1); 
202       pmSignal->SetMarkerColor(markerColorSignal);
203       pmSignal->SetMarkerStyle(1);
204       pmSignal->Draw();
205    
206     }
207
208   }
209
210   TGeometry *geoAlice = gAlice->GetGeometry();
211   TNode     *main     = (TNode *) ((geoAlice->GetListOfNodes())->First());
212   TIter      next(main->GetListOfNodes());
213   TNode     *module   = 0;
214   while ((module = (TNode *) next())) {
215     Char_t ch[100];
216     sprintf(ch,"%s\n",module->GetTitle());
217     if ((ch[0] == 'T') && ((ch[1] == 'R') || (ch[1] == 'P'))) {
218       module->SetVisibility( 3);
219     }
220     else {
221       module->SetVisibility(-1);
222     }
223   }
224       
225   geoAlice->Draw("same");
226
227   c1->Modified(); 
228   c1->Update(); 
229
230   printf("<AliTRDdisplayDigits3D> Number of digits:\n");
231   printf("                        signal = %d, bgnd = %d, merged = %d\n"
232         ,totalsignal,totalbgnd,totalmerged);
233
234   return 0;
235
236 }
237