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