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