]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RALICE/AliDevice.cxx
28-nov-2004 NvE User definable status word introduced in AliDevice to allow indicatio...
[u/mrichter/AliRoot.git] / RALICE / AliDevice.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 // $Id$
17
18 ///////////////////////////////////////////////////////////////////////////
19 // Class AliDevice
20 // Signal (Hit) handling of a generic device.
21 // Basically this class provides a user interface to group and handle
22 // various instances of AliSignal objects, called generically "hits".
23 // An AliDevice object itself has (in addition to hit storage) also the
24 // complete functionality of the class AliSignal.
25 //
26 // Example :
27 // =========
28 //
29 // AliDevice m;
30 // // Set user defined status word to indicate e.g. readout electronics version 
31 // m.SetStatus(100201);
32 // m.SetHitCopy(1);
33 // m.SetName("OM123");
34 //
35 // Float_t pos[3]={1,2,3};
36 // m.SetPosition(pos,"car");
37 //
38 // AliSignal s;
39 //
40 // s.Reset(1);
41 // s.SetName("OM123 Hit 1");
42 // s.SetSlotName("ADC");
43 // s.SetSignal(10);
44 // s.SetSlotName("LE",2);
45 // s.SetSignal(-100,2);
46 // s.SetSlotName("TOT",3);
47 // s.SetSignal(-1000,3);
48 // m.AddHit(s);
49 //
50 // s.Reset(1);
51 // s.SetName("OM123 Hit 2");
52 // s.SetSlotName("ADC");
53 // s.SetSignal(11);
54 // s.SetSlotName("LE",2);
55 // s.SetSignal(-101,2);
56 // s.SetSlotName("TOT",3);
57 // s.SetSignal(1001,3);
58 // m.AddHit(s);
59 //
60 // s.Reset(1);
61 // s.SetName("OM123 Hit 3");
62 // s.SetSlotName("ADC");
63 // s.SetSignal(12);
64 // s.SetSlotName("LE",2);
65 // s.SetSignal(-102,2);
66 // s.SetSlotName("TOT",3);
67 // s.SetSignal(-1002,3);
68 // m.AddHit(s);
69 //
70 // TObjArray* ordered=m.SortHits("TOT");
71 // nhits=ordered->GetEntries();
72 // for (Int_t i=0; i<nhits; i++)
73 // {
74 //  AliSignal* sx=(AliSignal*)ordered->At(i);
75 //  if (sx) sx->Data();
76 // }
77 //
78 //--- Author: Nick van Eijndhoven 23-jun-2004 Utrecht University
79 //- Modified: NvE $Date$ Utrecht University
80 ///////////////////////////////////////////////////////////////////////////
81
82 #include "AliDevice.h"
83 #include "Riostream.h"
84  
85 ClassImp(AliDevice) // Class implementation to enable ROOT I/O
86  
87 AliDevice::AliDevice() : AliSignal()
88 {
89 // Default constructor.
90 // The user definable status word is set to zero.
91 // By default private copies of the recorded hits will be made.
92 // This implies that by default the device will own the registered hits.
93 // See the SetHitCopy() memberfunction for further details.
94  fStatus=0;
95  fHitCopy=1;
96  fHits=0;
97  fOrdered=0;
98  fMarkers=0;
99 }
100 ///////////////////////////////////////////////////////////////////////////
101 AliDevice::~AliDevice()
102 {
103 // Default destructor.
104
105  // Remove backward links to this device from the hits
106  // which were not owned by it.
107  if (!fHitCopy)
108  {
109   for (Int_t ih=1; ih<=GetNhits(); ih++)
110   {
111    AliSignal* sx=GetHit(ih);
112    if (sx) sx->ResetLinks(this);
113   }
114  }
115
116  if (fHits)
117  {
118   delete fHits;
119   fHits=0;
120  }
121
122  if (fOrdered)
123  {
124   delete fOrdered;
125   fOrdered=0;
126  }
127
128  if (fMarkers)
129  {
130   delete fMarkers;
131   fMarkers=0;
132  }
133 }
134 ///////////////////////////////////////////////////////////////////////////
135 AliDevice::AliDevice(const AliDevice& dev) : AliSignal(dev)
136 {
137 // Copy constructor.
138
139  fHits=0;
140  fOrdered=0;
141  fMarkers=0;
142
143  fStatus=dev.GetStatus();
144  fHitCopy=dev.GetHitCopy();
145
146  Int_t nhits=dev.GetNhits();
147  if (nhits)
148  {
149   fHits=new TObjArray(nhits);
150   if (fHitCopy) fHits->SetOwner();
151   for (Int_t ih=1; ih<=nhits; ih++)
152   {
153    AliSignal* sx=dev.GetHit(ih);
154    if (fHitCopy)
155    {
156     fHits->Add(sx->Clone());
157     AliSignal* s=(AliSignal*)fHits->Last();
158     s->ResetLinks((AliDevice*)&dev);
159     s->SetDevice(this);
160    }
161    else
162    {
163     sx->AddLink(this);
164     fHits->Add(sx);
165    }
166   }
167  }
168 }
169 ///////////////////////////////////////////////////////////////////////////
170 void AliDevice::Reset(Int_t mode)
171 {
172 // Reset registered hits and AliSignal attributes.
173 // Note : The status word and HitCopy flag are NOT modified.
174 //        Use SetStatus() and SetHitCopy() to modify these parameters. 
175 // See AliSignal::Reset() for further details.
176  RemoveHits();
177  AliSignal::Reset(mode);
178 }
179 ///////////////////////////////////////////////////////////////////////////
180 void AliDevice::SetHitCopy(Int_t j)
181 {
182 // (De)activate the creation of private copies of the AliSignals added as hits.
183 // j=0 ==> No private copies are made; pointers of original hits are stored.
184 // j=1 ==> Private copies of the hits are made and these pointers are stored.
185 //
186 // Note : Once the storage contains pointer(s) to hit(s) one cannot
187 //        change the HitCopy mode anymore.
188 //        To change the HitCopy mode for an existing AliDevice containing
189 //        hits one first has to invoke either RemoveHits() or Reset().
190  if (!fHits)
191  {
192   if (j==0 || j==1)
193   {
194    fHitCopy=j;
195   }
196   else
197   {
198    cout << "*AliDevice::SetHitCopy* Invalid argument : " << j << endl;
199   }
200  }
201  else
202  {
203   cout << "*AliDevice::SetHitCopy* Storage already contained hits."
204        << "  ==> HitCopy mode not changed." << endl; 
205  }
206 }
207 ///////////////////////////////////////////////////////////////////////////
208 Int_t AliDevice::GetHitCopy() const
209 {
210 // Provide value of the HitCopy mode.
211 // 0 ==> No private copies are made; pointers of original hits are stored.
212 // 1 ==> Private copies of the hits are made and these pointers are stored.
213  return fHitCopy;
214 }
215 ///////////////////////////////////////////////////////////////////////////
216 void AliDevice::SetStatus(Int_t word)
217 {
218 // Set a user defined status word for this device.
219  fStatus=word;
220 }
221 ///////////////////////////////////////////////////////////////////////////
222 Int_t AliDevice::GetStatus() const
223 {
224 // Provide the user defined status word for this device.
225  return fStatus;
226 }
227 ///////////////////////////////////////////////////////////////////////////
228 void AliDevice::AddHit(AliSignal& s)
229 {
230 // Register an AliSignal object as a hit to this device.
231 // Note : In case this device owns the AliSignal object, the pointer to
232 //        this device will be stored in the special owning device
233 //        pointer of the AliSignal object.
234 //        In case this device does not own the AliSignal object, a (backward)
235 //        link to this device is added to the first slot of the AliSignal
236 //        if there was no link to this device already present.
237
238  if (!fHits)
239  {
240   fHits=new TObjArray(1);
241   if (fHitCopy) fHits->SetOwner();
242  }
243
244  // Check if this signal is already stored for this device.
245  Int_t nhits=GetNhits();
246  for (Int_t i=0; i<nhits; i++)
247  {
248   if (&s==fHits->At(i)) return; 
249  }
250
251  // Check for existing (backward) link to this device.
252  Int_t nlinks=s.GetNlinks(this);
253
254  if (fHitCopy)
255  {
256   fHits->Add(s.Clone());
257   // Remove unnecessary backward link(s) from the various slots
258   // and set the owning link to this device
259   AliSignal* sx=(AliSignal*)fHits->Last();
260   if (nlinks) sx->ResetLinks(this);
261   sx->SetDevice(this);
262  }
263  else
264  {
265   fHits->Add(&s);
266   // Set (backward) link to the this device
267   if (!nlinks) s.AddLink(this);
268  }
269 }
270 ///////////////////////////////////////////////////////////////////////////
271 void AliDevice::RemoveHit(AliSignal& s)
272 {
273 // Remove AliSignal object registered as a hit from this device.
274  if (fHits)
275  {
276   AliSignal* test=(AliSignal*)fHits->Remove(&s);
277   if (test)
278   {
279    fHits->Compress();
280    if (fHitCopy) delete test;
281   }
282  }
283  if (fOrdered)
284  {
285   AliSignal* test=(AliSignal*)fOrdered->Remove(&s);
286   if (test) fOrdered->Compress();
287  }
288 }
289 ///////////////////////////////////////////////////////////////////////////
290 void AliDevice::RemoveHits()
291 {
292 // Remove all AliSignal objects registered as hits from this device.
293  if (fHits)
294  {
295   delete fHits;
296   fHits=0;
297  }
298  if (fOrdered)
299  {
300   delete fOrdered;
301   fOrdered=0;
302  }
303  if (fMarkers)
304  {
305   delete fMarkers;
306   fMarkers=0;
307  }
308 }
309 ///////////////////////////////////////////////////////////////////////////
310 Int_t AliDevice::GetNhits() const
311 {
312 // Provide the number of registered hits for this device.
313  Int_t nhits=0;
314  if (fHits) nhits=fHits->GetEntries();
315  return nhits;
316 }
317 ///////////////////////////////////////////////////////////////////////////
318 AliSignal* AliDevice::GetHit(Int_t j) const
319 {
320 // Provide the AliSignal object registered as hit number j.
321 // Note : j=1 denotes the first hit.
322  if (!fHits) return 0;
323
324  if ((j >= 1) && (j <= GetNhits()))
325  {
326   return (AliSignal*)fHits->At(j-1);
327  }
328  else
329  {
330   return 0;
331  }
332 }
333 ///////////////////////////////////////////////////////////////////////////
334 AliSignal* AliDevice::GetIdHit(Int_t id) const
335 {
336 // Return the hit with unique identifier "id".
337  if (!fHits || id<0) return 0;
338
339  AliSignal* sx=0;
340  Int_t sid=0;
341  for (Int_t i=0; i<GetNhits(); i++)
342  {
343   sx=(AliSignal*)fHits->At(i);
344   if (sx)
345   {
346    sid=sx->GetUniqueID();
347    if (id==sid) return sx;
348   }
349  }
350  return 0; // No matching id found
351 }
352 ///////////////////////////////////////////////////////////////////////////
353 TObjArray* AliDevice::GetHits()
354 {
355 // Provide the references to all the registered hits.
356  return fHits;
357 }
358 ///////////////////////////////////////////////////////////////////////////
359 void AliDevice::ShowHit(Int_t j) const
360 {
361 // Show data of the registered j-th hit.
362 // If j=0 all associated hits will be shown.
363 // The default is j=0.
364  if (!j)
365  {
366   Int_t nhits=GetNhits();
367   for (Int_t ih=1; ih<=nhits; ih++)
368   {
369    AliSignal* sx=GetHit(ih);
370    if (sx) sx->Data();
371   }
372  }
373  else
374  {
375   AliSignal* s=GetHit(j);
376   if (s) s->Data();
377  }
378 }
379 ///////////////////////////////////////////////////////////////////////////
380 void AliDevice::Data(TString f) const
381 {
382 // Print the device and all registered hit info according to the specified
383 // coordinate frame.
384  AliSignal::Data(f);
385  Int_t nhits=GetNhits();
386  if (nhits)
387  {
388   cout << " The following " << nhits << " hits are registered : " << endl;
389   ShowHit();
390  }
391  else
392  {
393   cout << " No hits have been registered for this device." << endl;
394  }
395 }
396 ///////////////////////////////////////////////////////////////////////////
397 void AliDevice::GetExtremes(Float_t& vmin,Float_t& vmax,Int_t idx,TObjArray* hits) const
398 {
399 // Provide the min. and max. signal values of an array of hits.
400 // The input argument "idx" denotes the index of the signal slots to be investigated.
401 // The default is idx=1;
402 // In case hits=0 (default), the registered hits of the current device are used. 
403 // Signals which were declared as "Dead" will be rejected.
404 // The gain etc... corrected signals will be used in the process.
405
406  vmin=0;
407  vmax=0;
408
409  if (!hits) hits=fHits;
410  
411  if (idx<=0 || !hits) return;
412
413  Int_t nhits=hits->GetEntries();
414
415  Float_t sig=0;
416  for (Int_t i=0; i<nhits; i++)
417  {
418   AliSignal* sx=(AliSignal*)hits->At(i);
419
420   if (!sx) continue;
421   if (idx > sx->GetNvalues()) continue; // User specified slotindex out of range for this signal
422   if (sx->GetDeadValue(idx)) continue;  // Only take alive signals
423
424   sig=sx->GetSignal(idx,1);
425   if (i==0)
426   {
427    vmin=sig;
428    vmax=sig;
429   }
430   else
431   {
432    if (sig<vmin) vmin=sig;
433    if (sig>vmax) vmax=sig;
434   }
435  }
436 }
437 ///////////////////////////////////////////////////////////////////////////
438 void AliDevice::GetExtremes(Float_t& vmin,Float_t& vmax,TString name,TObjArray* hits) const
439 {
440 // Provide the min. and max. signal values of an array of hits.
441 // The input argument "name" denotes the name of the signal slots to be investigated.
442 // In case hits=0 (default), the registered hits of the current device are used. 
443 // Signals which were declared as "Dead" will be rejected.
444 // The gain etc... corrected signals will be used in the process.
445
446  vmin=0;
447  vmax=0;
448
449  if (!hits) hits=fHits;
450  
451  if (!hits) return;
452
453  Int_t nhits=hits->GetEntries();
454
455  Int_t idx=0; // The signal slotindex to perform the sorting on
456
457  Float_t sig=0;
458  for (Int_t i=0; i<nhits; i++)
459  {
460   AliSignal* sx=(AliSignal*)hits->At(i);
461
462   if (!sx) continue;
463
464   // Obtain the slotindex corresponding to the user selection
465   idx=sx->GetSlotIndex(name);
466   if (!idx) continue;
467
468   if (sx->GetDeadValue(idx)) continue; // Only take alive signals
469
470   sig=sx->GetSignal(idx,1);
471   if (i==0)
472   {
473    vmin=sig;
474    vmax=sig;
475   }
476   else
477   {
478    if (sig<vmin) vmin=sig;
479    if (sig>vmax) vmax=sig;
480   }
481  }
482 }
483 ///////////////////////////////////////////////////////////////////////////
484 TObjArray* AliDevice::SortHits(Int_t idx,Int_t mode,TObjArray* hits)
485 {
486 // Order the references to an array of hits by looping over the input array "hits"
487 // and checking the signal value. The ordered array is returned as a TObjArray.
488 // In case hits=0 (default), the registered hits of the current device are used. 
489 // Note that the original hit array is not modified.
490 // A "hit" represents an abstract object which is derived from AliSignal.
491 // The user can specify the index of the signal slot to perform the sorting on.
492 // By default the slotindex will be 1.
493 // Via the "mode" argument the user can specify ordering in decreasing
494 // order (mode=-1) or ordering in increasing order (mode=1).
495 // The default is mode=-1.
496 // Signals which were declared as "Dead" will be rejected.
497 // The gain etc... corrected signals will be used in the ordering process.
498
499  if (fOrdered)
500  {
501   delete fOrdered;
502   fOrdered=0;
503  }
504
505  if (!hits) hits=fHits;
506  
507  if (idx<=0 || abs(mode)!=1 || !hits) return fOrdered;
508
509  Int_t nhits=hits->GetEntries();
510  if (!nhits)
511  {
512   return fOrdered;
513  }
514  else
515  {
516   fOrdered=new TObjArray(nhits);
517  }
518
519  Int_t nord=0;
520  for (Int_t i=0; i<nhits; i++) // Loop over all hits of the array
521  {
522   AliSignal* s=(AliSignal*)hits->At(i);
523
524   if (!s) continue;
525
526   if (idx > s->GetNvalues()) continue; // User specified slotindex out of range for this signal
527   if (s->GetDeadValue(idx)) continue;  // Only take alive signals
528  
529   if (nord == 0) // store the first hit with a signal at the first ordered position
530   {
531    nord++;
532    fOrdered->AddAt(s,nord-1);
533    continue;
534   }
535  
536   for (Int_t j=0; j<=nord; j++) // put hit in the right ordered position
537   {
538    if (j == nord) // module has smallest (mode=-1) or largest (mode=1) signal seen so far
539    {
540     nord++;
541     fOrdered->AddAt(s,j); // add hit at the end
542     break; // go for next hit
543    }
544  
545    if (mode==-1 && s->GetSignal(idx,1) < ((AliSignal*)fOrdered->At(j))->GetSignal(idx,1)) continue;
546    if (mode==1 && s->GetSignal(idx,1) > ((AliSignal*)fOrdered->At(j))->GetSignal(idx,1)) continue;
547  
548    nord++;
549    for (Int_t k=nord-1; k>j; k--) // create empty position
550    {
551     fOrdered->AddAt(fOrdered->At(k-1),k);
552    }
553    fOrdered->AddAt(s,j); // put hit at empty position
554    break; // go for next matrix module
555   }
556  }
557  return fOrdered;
558 }
559 ///////////////////////////////////////////////////////////////////////////
560 TObjArray* AliDevice::SortHits(TString name,Int_t mode,TObjArray* hits)
561 {
562 // Order the references to an array of hits by looping over the input array "hits"
563 // and checking the signal value. The ordered array is returned as a TObjArray.
564 // In case hits=0 (default), the registered hits of the current device are used. 
565 // Note that the input array is not modified.
566 // A "hit" represents an abstract object which is derived from AliSignal.
567 // The user can specify the name of the signal slot to perform the sorting on.
568 // In case no matching slotname is found, the signal will be skipped.
569 // Via the "mode" argument the user can specify ordering in decreasing
570 // order (mode=-1) or ordering in increasing order (mode=1).
571 // The default is mode=-1.
572 // Signals which were declared as "Dead" will be rejected.
573 // The gain etc... corrected signals will be used in the ordering process.
574
575  if (fOrdered)
576  {
577   delete fOrdered;
578   fOrdered=0;
579  }
580
581  if (!hits) hits=fHits;
582  
583  if (abs(mode)!=1 || !hits) return fOrdered;
584
585  Int_t nhits=hits->GetEntries();
586  if (!nhits)
587  {
588   return fOrdered;
589  }
590  else
591  {
592   fOrdered=new TObjArray(nhits);
593  }
594
595  Int_t idx=0; // The signal slotindex to perform the sorting on
596
597  Int_t nord=0;
598  for (Int_t i=0; i<nhits; i++) // loop over all hits of the array
599  {
600   AliSignal* s=(AliSignal*)hits->At(i);
601
602   if (!s) continue;
603
604   // Obtain the slotindex corresponding to the user selection
605   idx=s->GetSlotIndex(name);
606   if (!idx) continue;
607
608   if (s->GetDeadValue(idx)) continue; // only take alive signals
609  
610   if (nord == 0) // store the first hit with a signal at the first ordered position
611   {
612    nord++;
613    fOrdered->AddAt(s,nord-1);
614    continue;
615   }
616  
617   for (Int_t j=0; j<=nord; j++) // put hit in the right ordered position
618   {
619    if (j == nord) // module has smallest (mode=-1) or largest (mode=1) signal seen so far
620    {
621     nord++;
622     fOrdered->AddAt(s,j); // add hit at the end
623     break; // go for next hit
624    }
625  
626    if (mode==-1 && s->GetSignal(idx,1) < ((AliSignal*)fOrdered->At(j))->GetSignal(idx,1)) continue;
627    if (mode==1 && s->GetSignal(idx,1) > ((AliSignal*)fOrdered->At(j))->GetSignal(idx,1)) continue;
628  
629    nord++;
630    for (Int_t k=nord-1; k>j; k--) // create empty position
631    {
632     fOrdered->AddAt(fOrdered->At(k-1),k);
633    }
634    fOrdered->AddAt(s,j); // put hit at empty position
635    break; // go for next matrix module
636   }
637  }
638  return fOrdered;
639 }
640 ///////////////////////////////////////////////////////////////////////////
641 void AliDevice::DisplayHits(Int_t idx,Float_t scale,TObjArray* hits,Int_t dp,Int_t mstyle,Int_t mcol)
642 {
643 // 3D color display of an array hits.
644 // The user can specify the index (default=1) of the signal slot to perform the display for.
645 // The marker size will indicate the absolute value of the signal (specified by the slotindex)
646 // as a percentage of the input argument "scale".
647 // In case scale<0 the maximum absolute signal value encountered in the hit array will be used
648 // to define the 100% scale. The default is scale=-1.
649 // In case hits=0 (default), the registered hits of the current device are used. 
650 // Note that the input array is not modified.
651 // In case dp=1 the device position will be used, otherwise the hit position will
652 // be used in the display. The default is dp=0.
653 // Via the "mstyle" and "mcol" arguments the user can specify the marker style
654 // and color (see TPolyMarker3D) respectively.
655 // The defaults are mstyle="large scalable dot" and mcol=blue.
656 // Signals which were declared as "Dead" will not be displayed.
657 // The gain etc... corrected signals will be used to determine the marker size.
658 //
659 // Note :
660 // ------
661 // Before any display activity, a TCanvas and a TView have to be initiated
662 // first by the user like for instance
663 // 
664 // TCanvas* c1=new TCanvas("c1","c1");
665 // TView* view=new TView(1);
666 // view->SetRange(-1000,-1000,-1000,1000,1000,1000);
667 // view->ShowAxis();
668
669  Int_t thisdev=0; // Indicate whether this is the owning device or not 
670  if (!hits)
671  {
672   hits=fHits;
673   thisdev=1;
674  }
675  
676  if (idx<=0 || !hits) return;
677
678  Int_t nhits=hits->GetEntries();
679  if (!nhits) return;
680
681  Float_t sigmax=fabs(scale);
682  if (scale<0)
683  {
684   Float_t vmin,vmax;
685   GetExtremes(vmin,vmax,idx,hits);
686   sigmax=fabs(vmax);
687   if (fabs(vmin)>sigmax) sigmax=fabs(vmin);
688  }
689
690  if (sigmax <=0) return;
691
692  if (fMarkers)
693  {
694   delete fMarkers;
695   fMarkers=0;
696  }
697  fMarkers=new TObjArray(nhits);
698  fMarkers->SetOwner();
699
700  Float_t pos[3];
701  GetPosition(pos,"car");
702
703  Float_t sig=0;
704  for (Int_t ih=0; ih<nhits; ih++)
705  {
706   AliSignal* sx=(AliSignal*)hits->At(ih);
707   if (!sx) continue;
708   if (!dp)
709   {
710    sx->GetPosition(pos,"car");
711   }
712   else
713   {
714    if (!thisdev)
715    {
716     AliDevice* dev=sx->GetDevice();
717     if (dev) dev->GetPosition(pos,"car");
718    }
719   }
720   sig=sx->GetSignal(idx,1);
721   TPolyMarker3D* m=new TPolyMarker3D();
722   m->SetMarkerStyle(mstyle);
723   m->SetMarkerColor(mcol);
724   m->SetMarkerSize(100.*fabs(sig)/sigmax);
725   m->SetPoint(0,pos[0],pos[1],pos[2]);
726   fMarkers->Add(m);
727   m->Draw();
728  }
729 }
730 ///////////////////////////////////////////////////////////////////////////
731 void AliDevice::DisplayHits(TString name,Float_t scale,TObjArray* hits,Int_t dp,Int_t mstyle,Int_t mcol)
732 {
733 // 3D color display of an array hits.
734 // The user can specify the name of the signal slot to perform the display for.
735 // The marker size will indicate the absolute value of the signal (specified by the slotname)
736 // as a percentage of the input argument "scale".
737 // In case scale<0 the maximum absolute signal value encountered in the hit array will be used
738 // to define the 100% scale. The default is scale=-1.
739 // In case hits=0 (default), the registered hits of the current device are used. 
740 // Note that the input array is not modified.
741 // In case dp=1 the device position will be used, otherwise the hit position will
742 // be used in the display. The default is dp=0.
743 // The marker size will indicate the percentage of the maximum encountered value
744 // of the absolute value of the name-specified input signal slots.
745 // Via the "mstyle" and "mcol" arguments the user can specify the marker style
746 // and color (see TPolyMarker3D) respectively.
747 // The defaults are mstyle="large scalable dot" and mcol=blue.
748 // Signals which were declared as "Dead" will not be displayed.
749 // The gain etc... corrected signals will be used to determine the marker size.
750 //
751 // Note :
752 // ------
753 // Before any display activity, a TCanvas and a TView have to be initiated
754 // first by the user like for instance
755 // 
756 // TCanvas* c1=new TCanvas("c1","c1");
757 // TView* view=new TView(1);
758 // view->SetRange(-1000,-1000,-1000,1000,1000,1000);
759 // view->ShowAxis();
760
761  Int_t thisdev=0; // Indicate whether this is the owning device or not 
762  if (!hits)
763  {
764   hits=fHits;
765   thisdev=1;
766  }
767  
768  if (!hits) return;
769
770  Int_t nhits=hits->GetEntries();
771
772  if (!nhits) return;
773
774  Float_t sigmax=fabs(scale);
775  if (scale<0)
776  {
777   Float_t vmin,vmax;
778   GetExtremes(vmin,vmax,name,hits);
779   sigmax=fabs(vmax);
780   if (fabs(vmin)>sigmax) sigmax=fabs(vmin);
781  }
782
783  if (sigmax <=0) return;
784
785  if (fMarkers)
786  {
787   delete fMarkers;
788   fMarkers=0;
789  }
790  fMarkers=new TObjArray(nhits);
791  fMarkers->SetOwner();
792
793  Float_t pos[3];
794  GetPosition(pos,"car");
795
796  Int_t idx=0; // The slot index corresponding to the user specified name
797  Float_t sig=0;
798  for (Int_t ih=0; ih<nhits; ih++)
799  {
800   AliSignal* sx=(AliSignal*)hits->At(ih);
801   if (!sx) continue;
802   idx=sx->GetSlotIndex(name);
803   if (!idx) continue;
804   if (!dp)
805   {
806    sx->GetPosition(pos,"car");
807   }
808   else
809   {
810    if (!thisdev)
811    {
812     AliDevice* dev=sx->GetDevice();
813     if (dev) dev->GetPosition(pos,"car");
814    }
815   }
816   sig=sx->GetSignal(idx,1);
817   TPolyMarker3D* m=new TPolyMarker3D();
818   m->SetMarkerStyle(mstyle);
819   m->SetMarkerColor(mcol);
820   m->SetMarkerSize(100.*fabs(sig)/sigmax);
821   m->SetPoint(0,pos[0],pos[1],pos[2]);
822   fMarkers->Add(m);
823   m->Draw();
824  }
825 }
826 ///////////////////////////////////////////////////////////////////////////
827 TObject* AliDevice::Clone(const char* name) const
828 {
829 // Make a deep copy of the current object and provide the pointer to the copy.
830 // This memberfunction enables automatic creation of new objects of the
831 // correct type depending on the object type, a feature which may be very useful
832 // for containers like AliEvent when adding objects in case the
833 // container owns the objects. This feature allows e.g. AliEvent
834 // to store either AliDevice objects or objects derived from AliDevice
835 // via tha AddDevice memberfunction, provided these derived classes also have
836 // a proper Clone memberfunction. 
837
838  AliDevice* dev=new AliDevice(*this);
839  if (name)
840  {
841   if (strlen(name)) dev->SetName(name);
842  }
843  return dev;
844 }
845 ///////////////////////////////////////////////////////////////////////////