]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RALICE/AliDevice.cxx
Updates on Lambda_c decays (S. Masciocchi)
[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::GetHit(TString name) const
335 {
336 // Provide the AliSignal object registered as hit with the specified name.
337 // Note : The first hit encountered with the specified name will be provided.
338
339  if (!fHits) return 0;
340
341  TString hitname;
342  Int_t nhits=GetNhits();
343  for (Int_t i=0; i<nhits; i++)
344  {
345   AliSignal* sx=(AliSignal*)fHits->At(i);
346   if (sx)
347   {
348    hitname=sx->GetName();
349    if (hitname == name) return sx;
350   }
351  }
352
353  return 0; // No matching name found
354 }
355 ///////////////////////////////////////////////////////////////////////////
356 AliSignal* AliDevice::GetIdHit(Int_t id) const
357 {
358 // Return the hit with unique identifier "id".
359  if (!fHits || id<0) return 0;
360
361  AliSignal* sx=0;
362  Int_t sid=0;
363  for (Int_t i=0; i<GetNhits(); i++)
364  {
365   sx=(AliSignal*)fHits->At(i);
366   if (sx)
367   {
368    sid=sx->GetUniqueID();
369    if (id==sid) return sx;
370   }
371  }
372  return 0; // No matching id found
373 }
374 ///////////////////////////////////////////////////////////////////////////
375 TObjArray* AliDevice::GetHits()
376 {
377 // Provide the references to all the registered hits.
378  return fHits;
379 }
380 ///////////////////////////////////////////////////////////////////////////
381 void AliDevice::ShowHit(Int_t j) const
382 {
383 // Show data of the registered j-th hit.
384 // If j=0 all associated hits will be shown.
385 // The default is j=0.
386  if (!j)
387  {
388   Int_t nhits=GetNhits();
389   for (Int_t ih=1; ih<=nhits; ih++)
390   {
391    AliSignal* sx=GetHit(ih);
392    if (sx) sx->Data();
393   }
394  }
395  else
396  {
397   AliSignal* s=GetHit(j);
398   if (s) s->Data();
399  }
400 }
401 ///////////////////////////////////////////////////////////////////////////
402 void AliDevice::Data(TString f,TString u) const
403 {
404 // Print the device and all registered hit info according to the specified
405 // coordinate frame f.
406 //
407 // The string argument "u" allows to choose between different angular units
408 // in case e.g. a spherical frame is selected.
409 // u = "rad" : angles provided in radians
410 //     "deg" : angles provided in degrees
411 //
412 // The defaults are f="car" and u="rad".
413
414  AliSignal::Data(f,u);
415  Int_t nhits=GetNhits();
416  if (nhits)
417  {
418   cout << " The following " << nhits << " hits are registered : " << endl;
419   ShowHit();
420  }
421  else
422  {
423   cout << " No hits have been registered for this device." << endl;
424  }
425 }
426 ///////////////////////////////////////////////////////////////////////////
427 void AliDevice::GetExtremes(Float_t& vmin,Float_t& vmax,Int_t idx,TObjArray* hits,Int_t mode) const
428 {
429 // Provide the min. and max. signal values of an array of hits.
430 // The input argument "idx" denotes the index of the signal slots to be investigated.
431 // The default is idx=1;
432 // In case hits=0 (default), the registered hits of the current device are used. 
433 // Signals which were declared as "Dead" will be rejected.
434 // The gain etc... corrected signals will be used in the process as specified
435 // by the  "mode" argument. The definition of this "mode" parameter corresponds to
436 // the description provided in the GetSignal memberfunction of class AliSignal.
437 // The default is mode=1 (for backward compatibility reasons).
438
439  vmin=0;
440  vmax=0;
441
442  if (!hits) hits=fHits;
443  
444  if (idx<=0 || !hits) return;
445
446  Int_t nhits=hits->GetEntries();
447
448  Float_t sig=0;
449  for (Int_t i=0; i<nhits; i++)
450  {
451   AliSignal* sx=(AliSignal*)hits->At(i);
452
453   if (!sx) continue;
454   if (idx > sx->GetNvalues()) continue; // User specified slotindex out of range for this signal
455   if (sx->GetDeadValue(idx)) continue;  // Only take alive signals
456
457   sig=sx->GetSignal(idx,mode);
458   if (i==0)
459   {
460    vmin=sig;
461    vmax=sig;
462   }
463   else
464   {
465    if (sig<vmin) vmin=sig;
466    if (sig>vmax) vmax=sig;
467   }
468  }
469 }
470 ///////////////////////////////////////////////////////////////////////////
471 void AliDevice::GetExtremes(Float_t& vmin,Float_t& vmax,TString name,TObjArray* hits,Int_t mode) const
472 {
473 // Provide the min. and max. signal values of an array of hits.
474 // The input argument "name" denotes the name of the signal slots to be investigated.
475 // In case hits=0 (default), the registered hits of the current device are used. 
476 // Signals which were declared as "Dead" will be rejected.
477 // The gain etc... corrected signals will be used in the process as specified
478 // by the  "mode" argument. The definition of this "mode" parameter corresponds to
479 // the description provided in the GetSignal memberfunction of class AliSignal.
480 // The default is mode=1 (for backward compatibility reasons).
481
482  vmin=0;
483  vmax=0;
484
485  if (!hits) hits=fHits;
486  
487  if (!hits) return;
488
489  Int_t nhits=hits->GetEntries();
490
491  Int_t idx=0; // The signal slotindex to perform the sorting on
492
493  Float_t sig=0;
494  for (Int_t i=0; i<nhits; i++)
495  {
496   AliSignal* sx=(AliSignal*)hits->At(i);
497
498   if (!sx) continue;
499
500   // Obtain the slotindex corresponding to the user selection
501   idx=sx->GetSlotIndex(name);
502   if (!idx) continue;
503
504   if (sx->GetDeadValue(idx)) continue; // Only take alive signals
505
506   sig=sx->GetSignal(idx,mode);
507   if (i==0)
508   {
509    vmin=sig;
510    vmax=sig;
511   }
512   else
513   {
514    if (sig<vmin) vmin=sig;
515    if (sig>vmax) vmax=sig;
516   }
517  }
518 }
519 ///////////////////////////////////////////////////////////////////////////
520 TObjArray* AliDevice::SortHits(Int_t idx,Int_t mode,TObjArray* hits,Int_t mcal)
521 {
522 // Order the references to an array of hits by looping over the input array "hits"
523 // and checking the signal value. The ordered array is returned as a TObjArray.
524 // In case hits=0 (default), the registered hits of the current device are used. 
525 // Note that the original hit array is not modified.
526 // A "hit" represents an abstract object which is derived from AliSignal.
527 // The user can specify the index of the signal slot to perform the sorting on.
528 // By default the slotindex will be 1.
529 // Via the "mode" argument the user can specify ordering in decreasing
530 // order (mode=-1) or ordering in increasing order (mode=1).
531 // The default is mode=-1.
532 // Signals which were declared as "Dead" will be rejected.
533 // The gain etc... corrected signals will be used in the ordering process as
534 // specified by the "mcal" argument. The definition of this "mcal" parameter
535 // corresponds to the signal correction mode described in the GetSignal
536 // memberfunction of class AliSignal.
537 // The default is mcal=1 (for backward compatibility reasons).
538
539  if (fOrdered)
540  {
541   delete fOrdered;
542   fOrdered=0;
543  }
544
545  if (!hits) hits=fHits;
546  
547  if (idx<=0 || abs(mode)!=1 || !hits) return fOrdered;
548
549  Int_t nhits=hits->GetEntries();
550  if (!nhits)
551  {
552   return fOrdered;
553  }
554  else
555  {
556   fOrdered=new TObjArray(nhits);
557  }
558
559  Int_t nord=0;
560  for (Int_t i=0; i<nhits; i++) // Loop over all hits of the array
561  {
562   AliSignal* s=(AliSignal*)hits->At(i);
563
564   if (!s) continue;
565
566   if (idx > s->GetNvalues()) continue; // User specified slotindex out of range for this signal
567   if (s->GetDeadValue(idx)) continue;  // Only take alive signals
568  
569   if (nord == 0) // store the first hit with a signal at the first ordered position
570   {
571    nord++;
572    fOrdered->AddAt(s,nord-1);
573    continue;
574   }
575  
576   for (Int_t j=0; j<=nord; j++) // put hit in the right ordered position
577   {
578    if (j == nord) // module has smallest (mode=-1) or largest (mode=1) signal seen so far
579    {
580     nord++;
581     fOrdered->AddAt(s,j); // add hit at the end
582     break; // go for next hit
583    }
584  
585    if (mode==-1 && s->GetSignal(idx,mcal) <= ((AliSignal*)fOrdered->At(j))->GetSignal(idx,mcal)) continue;
586    if (mode==1 && s->GetSignal(idx,mcal) >= ((AliSignal*)fOrdered->At(j))->GetSignal(idx,mcal)) continue;
587  
588    nord++;
589    for (Int_t k=nord-1; k>j; k--) // create empty position
590    {
591     fOrdered->AddAt(fOrdered->At(k-1),k);
592    }
593    fOrdered->AddAt(s,j); // put hit at empty position
594    break; // go for next hit
595   }
596  }
597  return fOrdered;
598 }
599 ///////////////////////////////////////////////////////////////////////////
600 TObjArray* AliDevice::SortHits(TString name,Int_t mode,TObjArray* hits,Int_t mcal)
601 {
602 // Order the references to an array of hits by looping over the input array "hits"
603 // and checking the signal value. The ordered array is returned as a TObjArray.
604 // In case hits=0 (default), the registered hits of the current device are used. 
605 // Note that the input array is not modified.
606 // A "hit" represents an abstract object which is derived from AliSignal.
607 // The user can specify the name of the signal slot to perform the sorting on.
608 // In case no matching slotname is found, the signal will be skipped.
609 // Via the "mode" argument the user can specify ordering in decreasing
610 // order (mode=-1) or ordering in increasing order (mode=1).
611 // The default is mode=-1.
612 // Signals which were declared as "Dead" will be rejected.
613 // The gain etc... corrected signals will be used in the ordering process as
614 // specified by the "mcal" argument. The definition of this "mcal" parameter
615 // corresponds to the signal correction mode described in the GetSignal
616 // memberfunction of class AliSignal.
617 // The default is mcal=1 (for backward compatibility reasons).
618
619  if (fOrdered)
620  {
621   delete fOrdered;
622   fOrdered=0;
623  }
624
625  if (!hits) hits=fHits;
626  
627  if (abs(mode)!=1 || !hits) return fOrdered;
628
629  Int_t nhits=hits->GetEntries();
630  if (!nhits)
631  {
632   return fOrdered;
633  }
634  else
635  {
636   fOrdered=new TObjArray(nhits);
637  }
638
639  Int_t idx=0; // The signal slotindex to perform the sorting on
640
641  Int_t nord=0;
642  for (Int_t i=0; i<nhits; i++) // loop over all hits of the array
643  {
644   AliSignal* s=(AliSignal*)hits->At(i);
645
646   if (!s) continue;
647
648   // Obtain the slotindex corresponding to the user selection
649   idx=s->GetSlotIndex(name);
650   if (!idx) continue;
651
652   if (s->GetDeadValue(idx)) continue; // only take alive signals
653  
654   if (nord == 0) // store the first hit with a signal at the first ordered position
655   {
656    nord++;
657    fOrdered->AddAt(s,nord-1);
658    continue;
659   }
660  
661   for (Int_t j=0; j<=nord; j++) // put hit in the right ordered position
662   {
663    if (j == nord) // module has smallest (mode=-1) or largest (mode=1) signal seen so far
664    {
665     nord++;
666     fOrdered->AddAt(s,j); // add hit at the end
667     break; // go for next hit
668    }
669  
670    if (mode==-1 && s->GetSignal(idx,mcal) <= ((AliSignal*)fOrdered->At(j))->GetSignal(idx,mcal)) continue;
671    if (mode==1 && s->GetSignal(idx,mcal) >= ((AliSignal*)fOrdered->At(j))->GetSignal(idx,mcal)) continue;
672  
673    nord++;
674    for (Int_t k=nord-1; k>j; k--) // create empty position
675    {
676     fOrdered->AddAt(fOrdered->At(k-1),k);
677    }
678    fOrdered->AddAt(s,j); // put hit at empty position
679    break; // go for next hit
680   }
681  }
682  return fOrdered;
683 }
684 ///////////////////////////////////////////////////////////////////////////
685 void AliDevice::DisplayHits(Int_t idx,Float_t scale,TObjArray* hits,Int_t dp,Int_t mode,Int_t mcol)
686 {
687 // 3D color display of an array hits.
688 // The user can specify the index (default=1) of the signal slot to perform the display for.
689 // The marker size will indicate the absolute value of the signal (specified by the slotindex)
690 // as a percentage of the input argument "scale".
691 // In case scale<0 the maximum absolute signal value encountered in the hit array will be used
692 // to define the 100% scale. The default is scale=-1.
693 // In case hits=0 (default), the registered hits of the current device are used. 
694 // Note that the input array is not modified.
695 // In case dp=1 the device position will be used, otherwise the hit position will
696 // be used in the display. The default is dp=0.
697 // Via the "mcol" argument the user can specify the marker color (see TPolyMarker3D).
698 // The default is mcol=blue.
699 // Signals which were declared as "Dead" will not be displayed.
700 // The gain etc... corrected signals will be used to determine the marker size.
701 // The gain correction is performed according to "mode" argument. The definition of this
702 // "mode" parameter corresponds to the description provided in the GetSignal
703 // memberfunction of class AliSignal.
704 // The default is mode=1 (for backward compatibility reasons).
705 //
706 // Note :
707 // ------
708 // Before any display activity, a TCanvas and a TView have to be initiated
709 // first by the user like for instance
710 // 
711 // TCanvas* c1=new TCanvas("c1","c1");
712 // TView* view=new TView(1);
713 // view->SetRange(-1000,-1000,-1000,1000,1000,1000);
714 // view->ShowAxis();
715
716  Int_t thisdev=0; // Indicate whether this is the owning device or not 
717  if (!hits)
718  {
719   hits=fHits;
720   thisdev=1;
721  }
722  
723  if (idx<=0 || !hits) return;
724
725  Int_t nhits=hits->GetEntries();
726  if (!nhits) return;
727
728  Float_t sigmax=fabs(scale);
729  if (scale<0)
730  {
731   Float_t vmin,vmax;
732   GetExtremes(vmin,vmax,idx,hits,mode);
733   sigmax=fabs(vmax);
734   if (fabs(vmin)>sigmax) sigmax=fabs(vmin);
735  }
736
737  if (sigmax <=0) return;
738
739  if (fMarkers)
740  {
741   delete fMarkers;
742   fMarkers=0;
743  }
744  fMarkers=new TObjArray(nhits);
745  fMarkers->SetOwner();
746
747  Float_t pos[3];
748  GetPosition(pos,"car");
749
750  Float_t sig=0;
751  for (Int_t ih=0; ih<nhits; ih++)
752  {
753   AliSignal* sx=(AliSignal*)hits->At(ih);
754   if (!sx) continue;
755   if (!dp)
756   {
757    sx->GetPosition(pos,"car");
758   }
759   else
760   {
761    if (!thisdev)
762    {
763     AliDevice* dev=sx->GetDevice();
764     if (dev) dev->GetPosition(pos,"car");
765    }
766   }
767   sig=sx->GetSignal(idx,mode);
768
769   // Skip dead signals
770   if (fabs(sig) <= 0.) continue;
771
772   TPolyMarker3D* m=new TPolyMarker3D();
773   m->SetMarkerStyle(8);
774   m->SetMarkerColor(mcol);
775   m->SetMarkerSize(100.*fabs(sig)/sigmax);
776   m->SetPoint(0,pos[0],pos[1],pos[2]);
777   fMarkers->Add(m);
778   m->Draw();
779  }
780 }
781 ///////////////////////////////////////////////////////////////////////////
782 void AliDevice::DisplayHits(TString name,Float_t scale,TObjArray* hits,Int_t dp,Int_t mode,Int_t mcol)
783 {
784 // 3D color display of an array hits.
785 // The user can specify the name of the signal slot to perform the display for.
786 // The marker size will indicate the absolute value of the signal (specified by the slotname)
787 // as a percentage of the input argument "scale".
788 // In case scale<0 the maximum absolute signal value encountered in the hit array will be used
789 // to define the 100% scale. The default is scale=-1.
790 // In case hits=0 (default), the registered hits of the current device are used. 
791 // Note that the input array is not modified.
792 // In case dp=1 the device position will be used, otherwise the hit position will
793 // be used in the display. The default is dp=0.
794 // The marker size will indicate the percentage of the maximum encountered value
795 // of the absolute value of the name-specified input signal slots.
796 // Via the "mcol" argument the user can specify the marker color (see TPolyMarker3D).
797 // The default is mcol=blue.
798 // Signals which were declared as "Dead" will not be displayed.
799 // The gain etc... corrected signals will be used to determine the marker size.
800 // The gain correction is performed according to "mode" argument. The definition of this
801 // "mode" parameter corresponds to the description provided in the GetSignal
802 // memberfunction of class AliSignal.
803 // The default is mode=1 (for backward compatibility reasons).
804 //
805 // Note :
806 // ------
807 // Before any display activity, a TCanvas and a TView have to be initiated
808 // first by the user like for instance
809 // 
810 // TCanvas* c1=new TCanvas("c1","c1");
811 // TView* view=new TView(1);
812 // view->SetRange(-1000,-1000,-1000,1000,1000,1000);
813 // view->ShowAxis();
814
815  Int_t thisdev=0; // Indicate whether this is the owning device or not 
816  if (!hits)
817  {
818   hits=fHits;
819   thisdev=1;
820  }
821  
822  if (!hits) return;
823
824  Int_t nhits=hits->GetEntries();
825
826  if (!nhits) return;
827
828  Float_t sigmax=fabs(scale);
829  if (scale<0)
830  {
831   Float_t vmin,vmax;
832   GetExtremes(vmin,vmax,name,hits,mode);
833   sigmax=fabs(vmax);
834   if (fabs(vmin)>sigmax) sigmax=fabs(vmin);
835  }
836
837  if (sigmax <=0) return;
838
839  if (fMarkers)
840  {
841   delete fMarkers;
842   fMarkers=0;
843  }
844  fMarkers=new TObjArray(nhits);
845  fMarkers->SetOwner();
846
847  Float_t pos[3];
848  GetPosition(pos,"car");
849
850  Int_t idx=0; // The slot index corresponding to the user specified name
851  Float_t sig=0;
852  for (Int_t ih=0; ih<nhits; ih++)
853  {
854   AliSignal* sx=(AliSignal*)hits->At(ih);
855   if (!sx) continue;
856   idx=sx->GetSlotIndex(name);
857   if (!idx) continue;
858   if (!dp)
859   {
860    sx->GetPosition(pos,"car");
861   }
862   else
863   {
864    if (!thisdev)
865    {
866     AliDevice* dev=sx->GetDevice();
867     if (dev) dev->GetPosition(pos,"car");
868    }
869   }
870   sig=sx->GetSignal(idx,mode);
871
872   // Skip dead signals
873   if (fabs(sig) <= 0.) continue;
874
875   TPolyMarker3D* m=new TPolyMarker3D();
876   m->SetMarkerStyle(8);
877   m->SetMarkerColor(mcol);
878   m->SetMarkerSize(100.*fabs(sig)/sigmax);
879   m->SetPoint(0,pos[0],pos[1],pos[2]);
880   fMarkers->Add(m);
881   m->Draw();
882  }
883 }
884 ///////////////////////////////////////////////////////////////////////////
885 TObject* AliDevice::Clone(const char* name) const
886 {
887 // Make a deep copy of the current object and provide the pointer to the copy.
888 // This memberfunction enables automatic creation of new objects of the
889 // correct type depending on the object type, a feature which may be very useful
890 // for containers like AliEvent when adding objects in case the
891 // container owns the objects. This feature allows e.g. AliEvent
892 // to store either AliDevice objects or objects derived from AliDevice
893 // via tha AddDevice memberfunction, provided these derived classes also have
894 // a proper Clone memberfunction. 
895
896  AliDevice* dev=new AliDevice(*this);
897  if (name)
898  {
899   if (strlen(name)) dev->SetName(name);
900  }
901  return dev;
902 }
903 ///////////////////////////////////////////////////////////////////////////