]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/Alieve/TRDModuleImp.cxx
Runloader is updated when moving to next file (quick fix).
[u/mrichter/AliRoot.git] / EVE / Alieve / TRDModuleImp.cxx
1 #include "TRDModuleImp.h"
2 #include "TRDData.h"
3
4 #include "TMath.h"
5 #include <TGListTree.h>
6
7 #include "Reve/RGTopFrame.h"
8 #include "Reve/Track.h"
9
10 #include "AliRun.h"
11 #include "AliTRDv1.h"
12 #include "AliTRDgeometry.h"
13 #include "AliTRDCommonParam.h"
14 #include "AliTRDpadPlane.h"
15 #include "AliTRDdigit.h"
16 #include "AliTRDhit.h"
17 #include "AliTRDcluster.h"
18 #include "AliTRDcalibDB.h"
19 #include "AliTRDdataArrayI.h"
20 #include "AliTRDmcmTracklet.h"
21
22
23
24 using namespace Reve;
25 using namespace Alieve;
26 using namespace std;
27
28 ClassImp(TRDChamber)
29 ClassImp(TRDNode)
30
31 //________________________________________________________
32 TRDNode::TRDNode(const char *typ, Int_t det) :
33   Reve::RenderElement(), TRDModule(typ, det)
34 {
35 }
36
37 //________________________________________________________
38 void    TRDNode::Paint(Option_t* option)
39 {
40         List_i iter = fChildren.begin();
41         while(iter != fChildren.end()){
42                 (dynamic_cast<TRDModule*>(*iter))->Paint(option);
43                 iter++;
44         }
45 }
46
47 //________________________________________________________
48 void    TRDNode::Reset()
49 {
50         List_i iter = fChildren.begin();
51         while(iter != fChildren.end()){
52                 (dynamic_cast<TRDModule*>(*iter))->Reset();
53                 iter++;
54         }
55 }
56
57 //________________________________________________________
58 void TRDNode::Collapse()
59 {
60         TGListTree *list = gReve->GetListTree();
61         TRDNode *node = 0x0;
62         List_i iter = fChildren.begin();
63         while(iter != fChildren.end()){
64                 if((node = dynamic_cast<TRDNode*>(*iter))) node->Collapse();
65                 list->CloseItem(FindListTreeItem(list));
66                 iter++;
67         }
68 }
69
70 //________________________________________________________
71 void TRDNode::Expand()
72 {
73         TGListTree *list = gReve->GetListTree();
74         TRDNode *node = 0x0;
75         List_i iter = fChildren.begin();
76         while(iter != fChildren.end()){
77                 if((node = dynamic_cast<TRDNode*>(*iter))) node->Expand();
78                 list->OpenItem(FindListTreeItem(list));
79                 iter++;
80         }
81 }
82
83 //________________________________________________________
84 void TRDNode::EnableListElements()
85 {
86         SetRnrSelf(kTRUE);
87         TRDNode *node = 0x0;
88         TRDChamber *chmb = 0x0; 
89         List_i iter = fChildren.begin();
90         while(iter != fChildren.end()){
91                 if((node = dynamic_cast<TRDNode*>(*iter))){
92                         node->SetRnrSelf(kTRUE);
93                         node->EnableListElements();
94                 }
95                 if((chmb = dynamic_cast<TRDChamber*>(*iter))) chmb->SetRnrSelf(kTRUE);
96                 iter++;
97         }
98         gReve->Redraw3D();
99 }
100
101 //________________________________________________________
102 void TRDNode::DisableListElements()
103 {
104         SetRnrSelf(kFALSE);
105         TRDNode *node = 0x0;
106         TRDChamber *chmb = 0x0; 
107         List_i iter = fChildren.begin();
108         while(iter != fChildren.end()){
109                 if((node = dynamic_cast<TRDNode*>(*iter))){
110                         node->SetRnrSelf(kFALSE);
111                         node->DisableListElements();
112                 }
113                 if((chmb = dynamic_cast<TRDChamber*>(*iter))) chmb->SetRnrSelf(kFALSE);
114                 iter++;
115         }
116         gReve->Redraw3D();
117 }
118
119 //________________________________________________________
120 void TRDNode::UpdateLeaves()
121 {
122         TRDModule *module;
123         List_i iter = fChildren.begin();
124         while(iter != fChildren.end()){
125                 module = dynamic_cast<TRDModule*>(*iter);
126                 if(!module) continue;
127                 
128                 module->fRnrHits = fRnrHits;
129                 module->fRnrDigits = fRnrDigits;
130                 module->fDigitsLog = fDigitsLog;
131                 module->fDigitsBox = fDigitsBox;
132                 module->fDigitsThreshold = fDigitsThreshold;
133                 module->kDigitsNeedRecompute = kDigitsNeedRecompute;
134                 module->fRnrRecPoints = fRnrRecPoints;
135                 module->fRnrTracklets = fRnrTracklets;
136                 iter++;
137         }
138
139         TRDNode *node = 0x0;
140         iter = fChildren.begin();
141         while(iter != fChildren.end()){
142                 if((node = dynamic_cast<TRDNode*>(*iter))) node->UpdateLeaves();
143                 iter++;
144         }
145 }
146
147
148 //________________________________________________________
149 void TRDNode::UpdateNode()
150 {
151 //      Info("UpdateNode()", Form("%s", GetName()));
152         TRDNode *node = 0x0;
153         List_i iter = fChildren.begin();
154         while(iter != fChildren.end()){
155                 if((node = dynamic_cast<TRDNode*>(*iter))) node->UpdateNode();
156                 iter++;
157         }
158
159         Int_t score[11];
160         for(int i=0; i<11; i++) score[i] = 0;
161         TRDModule *module;
162         iter = fChildren.begin();
163         while(iter != fChildren.end()){
164                 module = dynamic_cast<TRDModule*>(*iter);
165                 if(!module) continue;
166                 score[0] += (module->fLoadHits) ? 1 : 0;
167                 score[1] += (module->fRnrHits) ? 1 : 0;
168                 
169                 score[2] += (module->fLoadDigits) ? 1 : 0;
170                 score[3] += (module->fRnrDigits) ? 1 : 0;
171                 score[4] += (module->fDigitsLog) ? 1 : 0;
172                 score[5] += (module->fDigitsBox) ? 1 : 0;
173                 score[6] += (module->kDigitsNeedRecompute) ? 1 : 0;
174
175                 score[7] += (module->fLoadRecPoints) ? 1 : 0;
176                 score[8] += (module->fRnrRecPoints) ? 1 : 0;
177
178                 score[9] += (module->fLoadTracklets) ? 1 : 0;
179                 score[10] += (module->fRnrTracklets) ? 1 : 0;
180                 iter++;
181         }
182
183         Int_t size = fChildren.size(); 
184         fLoadHits      = (score[0] > 0) ? kTRUE : kFALSE;
185         fRnrHits       = (score[1] == size) ? kTRUE : kFALSE;
186
187         fLoadDigits    = (score[2] > 0) ? kTRUE : kFALSE;
188         fRnrDigits     = (score[3] == size) ? kTRUE : kFALSE;
189         fDigitsLog     = (score[4] == size) ? kTRUE : kFALSE;
190         fDigitsBox     = (score[5] == size) ? kTRUE : kFALSE;
191         kDigitsNeedRecompute = (score[6] == size) ? kTRUE : kFALSE;
192
193         fLoadRecPoints = (score[7] > 0) ? kTRUE : kFALSE;
194         fRnrRecPoints  = (score[8] == size) ? kTRUE : kFALSE;
195
196         fLoadTracklets = (score[9] > 0) ? kTRUE : kFALSE;
197         fRnrTracklets  = (score[10] == size) ? kTRUE : kFALSE;
198 }
199
200
201 ///////////////////////////////////////////////////////////
202 /////////////        TRDChamber       /////////////////////
203 ///////////////////////////////////////////////////////////
204
205 //________________________________________________________
206 TRDChamber::TRDChamber(Int_t det) :
207   Reve::RenderElement(), TRDModule("Chmb", det), rowMax(-1), colMax(-1), timeMax(22), fX0(0.), fPla(-1)
208 {
209   //
210   // Constructor
211   //
212         
213         fDigits    = 0x0;
214         fHits      = 0x0;
215         fRecPoints = 0x0;
216         fTracklets = 0x0;
217         
218         AliTRDCommonParam* parCom = AliTRDCommonParam::Instance();
219         samplingFrequency = parCom->GetSamplingFrequency();
220         
221         fGeo      = 0x0;
222         fPadPlane = 0x0;
223 }
224
225 //________________________________________________________
226 TRDChamber::TRDChamber(const TRDChamber &mod):
227   Reve::RenderElement(), TRDModule("Chmb", mod.fDet)
228 {
229   //
230   // Copy constructor
231   //
232         SetMainColor(mod.GetMainColor());
233
234         if(mod.fDigits) {}
235         if(mod.fHits) {}
236         if(mod.fRecPoints){} 
237 }
238
239 //________________________________________________________
240 TRDChamber& TRDChamber::operator=(const TRDChamber &mod)
241 {
242   //
243   // Assignment operator
244   //
245
246   if (this != &mod) {
247     fDet    = mod.fDet;
248                 if(mod.fDigits) {}
249                 if(mod.fHits) {}
250                 if(mod.fRecPoints){} 
251   }
252   return *this;
253 }
254
255 //________________________________________________________
256 void TRDChamber::LoadClusters(TObjArray *clusters)
257 {
258   //
259   // Draw clusters
260   //
261         
262         if(!fGeo){
263                 Error("LoadClusters()", Form("Geometry not set for chamber %d. Please call first TRDChamber::SetGeometry().", fDet));
264                 return;
265         }
266 //      Info("LoadClusters()", Form("clusters = 0x%x", clusters));
267         if(!fRecPoints){
268                 fRecPoints = new TRDHits("clusters");
269                 fRecPoints->SetMarkerSize(1.);
270                 fRecPoints->SetMarkerStyle(24);
271                 fRecPoints->SetMarkerColor(6);
272         } else fRecPoints->Reset();
273
274         Float_t q, z0;
275   Double_t cloc[3], cglo[3];
276         
277         z0 = fGeo->GetTime0(fPla);      
278         AliTRDcluster *c=0x0;
279         for(int iclus=0; iclus<clusters->GetEntriesFast(); iclus++){
280                 c = (AliTRDcluster*)clusters->UncheckedAt(iclus);
281                 cloc[2] = c->GetZ(); //x
282                 cloc[1] = c->GetY(); //y
283                 cloc[0] = z0 - c->GetX(); //z
284                 q = c->GetQ();
285                 fGeo->RotateBack(fDet,cloc,cglo);
286                 fRecPoints->SetNextPoint(cglo[0], cglo[1], cglo[2]);
287                 fRecPoints->SetPointId(this);
288         }
289         fLoadRecPoints = kTRUE;
290 }
291
292 //________________________________________________________
293 void TRDChamber::LoadDigits(AliTRDdigitsManager *digits)
294 {
295   //
296   // Draw digits
297   //
298         if(!fGeo){
299                 Error("LoadDigits()", Form("Geometry not set for chamber %d. Please call first TRDChamber::SetGeometry().", fDet));
300                 return;
301         }
302 //      Info("LoadDigits()", Form("digits =0x%x", digits));
303         
304         if(!fDigits) fDigits = new TRDDigits(this);
305         else fDigits->Reset();
306         
307         fDigits->SetData(digits);
308         fLoadDigits = kTRUE;
309 }
310
311 //________________________________________________________
312 void TRDChamber::AddHit(AliTRDhit *hit)
313 {
314   //
315   // Draw hits
316   //
317 //      Info("AddHit()", Form("%s", GetName()));
318
319         if(!fHits){
320                 fHits = new TRDHits("hits");
321                 fHits->SetMarkerSize(.1);
322                 fHits->SetMarkerColor(2);
323         }
324         
325         fHits->SetNextPoint(hit->X(), hit->Y(), hit->Z());
326         fHits->SetPointId(this);
327         fLoadHits = kTRUE;
328 }
329
330 //________________________________________________________
331 void TRDChamber::LoadTracklets(TObjArray *tracks)
332 {
333   //
334   // Draw tracks
335   //
336         if(!fGeo){
337                 Error("LoadTracklets()", Form("Geometry not set for chamber %d. Please call first TRDChamber::SetGeometry().", fDet));
338                 return;
339         }
340 //      Info("LoadTracklets()", Form("tracks = 0x%x", tracks));
341         
342         if(!fTracklets){
343                 fTracklets = new std::vector<Reve::Track*>;
344         } else fTracklets->clear();
345         
346         
347         AliTRDmcmTracklet *trk = 0x0;
348         Double_t cloc[3], cglo[3];
349         for(int itrk=0; itrk<tracks->GetEntries();itrk++){
350                 trk = (AliTRDmcmTracklet*)tracks->At(itrk);
351                 trk->MakeTrackletGraph(fGeo,.5);
352                 fTracklets->push_back(new Reve::Track());
353                 fTracklets->back()->SetLineColor(4);
354                 
355                 cloc[0] = trk->GetTime0(); // x0
356                 cloc[1] = trk->GetOffset(); // y0
357                 cloc[2] = trk->GetRowz(); // z
358           fGeo->RotateBack(fDet,cloc,cglo);
359                 fTracklets->back()->SetNextPoint(cglo[0], cglo[1], cglo[2]);
360                 
361                 cloc[0] += 3.7; // x1
362                 cloc[1] += TMath::Tan(trk->GetSlope()*TMath::Pi()/180.) * 3.7; // y1
363           fGeo->RotateBack(fDet,cloc,cglo);
364                 fTracklets->back()->SetNextPoint(cglo[0], cglo[1], cglo[2]);
365         }
366         fLoadTracklets = kTRUE;
367 }
368
369 //____________________________________________________
370 void    TRDChamber::Paint(Option_t* option)
371 {
372 /*      Info("Paint()", Form("%s", GetName()));*/
373         if(!fRnrSelf) return;
374         if(fDigits && fRnrDigits){
375                 if(kDigitsNeedRecompute){
376                         fDigits->ComputeRepresentation();
377                         kDigitsNeedRecompute = kFALSE;
378                 }
379                 fDigits->Paint(option);
380         }
381         if(fRecPoints && fRnrRecPoints) fRecPoints->GetObject()->Paint(option);
382         if(fHits && fRnrHits) fHits->GetObject()->Paint(option);
383         if(fTracklets && fRnrTracklets){
384                 for(vector<Reve::Track*>::iterator i=fTracklets->begin(); i != fTracklets->end(); ++i) (*i)->Paint(option);
385         }
386 }
387
388 //________________________________________________________
389 void    TRDChamber::Reset()
390 {
391         if(fHits){
392                 fHits->Reset();
393                 fLoadHits = kFALSE;
394         }
395         if(fDigits){
396                 fDigits->Reset();
397                 fLoadDigits = kFALSE;
398         }
399         if(fRecPoints){
400                 fRecPoints->Reset();
401                 fLoadRecPoints = kFALSE;
402         }
403         if(fTracklets){
404                 fTracklets->clear();
405                 fLoadTracklets = kFALSE;
406         }
407 }
408
409 //________________________________________________________
410 void TRDChamber::SetGeometry(AliTRDgeometry *geo)
411 {
412         fGeo = geo;
413                 
414         fPla = fGeo->GetPlane(fDet);
415         fX0 = fGeo->GetTime0(fPla);
416         
417         AliTRDCommonParam *parcom = AliTRDCommonParam::Instance();
418         fPadPlane = parcom->GetPadPlane(fPla,fGeo->GetChamber(fDet));
419         rowMax = fPadPlane->GetNrows();
420         colMax = fPadPlane->GetNcols();
421 }
422