1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////
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.
31 // m.SetName("OM123");
33 // Float_t pos[3]={1,2,3};
34 // m.SetPosition(pos,"car");
39 // s.SetName("OM123 Hit 1");
40 // s.SetSlotName("ADC");
42 // s.SetSlotName("LE",2);
43 // s.SetSignal(-100,2);
44 // s.SetSlotName("TOT",3);
45 // s.SetSignal(-1000,3);
49 // s.SetName("OM123 Hit 2");
50 // s.SetSlotName("ADC");
52 // s.SetSlotName("LE",2);
53 // s.SetSignal(-101,2);
54 // s.SetSlotName("TOT",3);
55 // s.SetSignal(1001,3);
59 // s.SetName("OM123 Hit 3");
60 // s.SetSlotName("ADC");
62 // s.SetSlotName("LE",2);
63 // s.SetSignal(-102,2);
64 // s.SetSlotName("TOT",3);
65 // s.SetSignal(-1002,3);
68 // TObjArray* ordered=m.SortHits("TOT");
69 // nhits=ordered->GetEntries();
70 // for (Int_t i=0; i<nhits; i++)
72 // AliSignal* sx=(AliSignal*)ordered->At(i);
73 // if (sx) sx->Data();
76 //--- Author: Nick van Eijndhoven 23-jun-2004 Utrecht University
77 //- Modified: NvE $Date$ Utrecht University
78 ///////////////////////////////////////////////////////////////////////////
80 #include "AliDevice.h"
81 #include "Riostream.h"
83 ClassImp(AliDevice) // Class implementation to enable ROOT I/O
85 AliDevice::AliDevice() : AliSignal()
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.
96 ///////////////////////////////////////////////////////////////////////////
97 AliDevice::~AliDevice()
99 // Default destructor.
101 // Remove backward links to this device from the hits
102 // which were not owned by it.
105 for (Int_t ih=1; ih<=GetNhits(); ih++)
107 AliSignal* sx=GetHit(ih);
108 if (sx) sx->ResetLinks(this);
130 ///////////////////////////////////////////////////////////////////////////
131 AliDevice::AliDevice(const AliDevice& dev) : AliSignal(dev)
139 fHitCopy=dev.GetHitCopy();
141 Int_t nhits=dev.GetNhits();
144 fHits=new TObjArray(nhits);
145 if (fHitCopy) fHits->SetOwner();
146 for (Int_t ih=1; ih<=nhits; ih++)
148 AliSignal* sx=dev.GetHit(ih);
151 fHits->Add(sx->Clone());
152 AliSignal* s=(AliSignal*)fHits->Last();
153 s->ResetLinks((AliDevice*)&dev);
164 ///////////////////////////////////////////////////////////////////////////
165 void AliDevice::Reset(Int_t mode)
167 // Reset registered hits and AliSignal attributes.
168 // See AliSignal::Reset() for further details.
170 AliSignal::Reset(mode);
172 ///////////////////////////////////////////////////////////////////////////
173 void AliDevice::SetHitCopy(Int_t j)
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.
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().
191 cout << "*AliDevice::SetHitCopy* Invalid argument : " << j << endl;
196 cout << "*AliDevice::SetHitCopy* Storage already contained hits."
197 << " ==> HitCopy mode not changed." << endl;
200 ///////////////////////////////////////////////////////////////////////////
201 Int_t AliDevice::GetHitCopy() const
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.
208 ///////////////////////////////////////////////////////////////////////////
209 void AliDevice::AddHit(AliSignal& s)
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.
221 fHits=new TObjArray(1);
222 if (fHitCopy) fHits->SetOwner();
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++)
229 if (&s==fHits->At(i)) return;
232 // Check for existing (backward) link to this device.
233 Int_t nlinks=s.GetNlinks(this);
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);
247 // Set (backward) link to the this device
248 if (!nlinks) s.AddLink(this);
251 ///////////////////////////////////////////////////////////////////////////
252 void AliDevice::RemoveHit(AliSignal& s)
254 // Remove AliSignal object registered as a hit from this device.
257 AliSignal* test=(AliSignal*)fHits->Remove(&s);
261 if (fHitCopy) delete test;
266 AliSignal* test=(AliSignal*)fOrdered->Remove(&s);
267 if (test) fOrdered->Compress();
270 ///////////////////////////////////////////////////////////////////////////
271 void AliDevice::RemoveHits()
273 // Remove all AliSignal objects registered as hits from this device.
290 ///////////////////////////////////////////////////////////////////////////
291 Int_t AliDevice::GetNhits() const
293 // Provide the number of registered hits for this device.
295 if (fHits) nhits=fHits->GetEntries();
298 ///////////////////////////////////////////////////////////////////////////
299 AliSignal* AliDevice::GetHit(Int_t j) const
301 // Provide the AliSignal object registered as hit number j.
302 // Note : j=1 denotes the first hit.
303 if (!fHits) return 0;
305 if ((j >= 1) && (j <= GetNhits()))
307 return (AliSignal*)fHits->At(j-1);
314 ///////////////////////////////////////////////////////////////////////////
315 TObjArray* AliDevice::GetHits()
317 // Provide the references to all the registered hits.
320 ///////////////////////////////////////////////////////////////////////////
321 void AliDevice::ShowHit(Int_t j) const
323 // Show data of the registered j-th hit.
324 // If j=0 all associated hits will be shown.
325 // The default is j=0.
328 Int_t nhits=GetNhits();
329 for (Int_t ih=1; ih<=nhits; ih++)
331 AliSignal* sx=GetHit(ih);
337 AliSignal* s=GetHit(j);
341 ///////////////////////////////////////////////////////////////////////////
342 void AliDevice::Data(TString f) const
344 // Print the device and all registered hit info according to the specified
347 Int_t nhits=GetNhits();
350 cout << " The following " << nhits << " hits are registered : " << endl;
355 cout << " No hits have been registered for this device." << endl;
358 ///////////////////////////////////////////////////////////////////////////
359 void AliDevice::GetExtremes(Float_t& vmin,Float_t& vmax,Int_t idx,TObjArray* hits) const
361 // Provide the min. and max. signal values of an array of hits.
362 // The input argument "idx" denotes the index of the signal slots to be investigated.
363 // The default is idx=1;
364 // In case hits=0 (default), the registered hits of the current device are used.
365 // Signals which were declared as "Dead" will be rejected.
366 // The gain etc... corrected signals will be used in the process.
371 if (!hits) hits=fHits;
373 if (idx<=0 || !hits) return;
375 Int_t nhits=hits->GetEntries();
378 for (Int_t i=0; i<nhits; i++)
380 AliSignal* sx=(AliSignal*)hits->At(i);
383 if (idx > sx->GetNvalues()) continue; // User specified slotindex out of range for this signal
384 if (sx->GetDeadValue(idx)) continue; // Only take alive signals
386 sig=sx->GetSignal(idx,1);
394 if (sig<vmin) vmin=sig;
395 if (sig>vmax) vmax=sig;
399 ///////////////////////////////////////////////////////////////////////////
400 void AliDevice::GetExtremes(Float_t& vmin,Float_t& vmax,TString name,TObjArray* hits) const
402 // Provide the min. and max. signal values of an array of hits.
403 // The input argument "name" denotes the name of the signal slots to be investigated.
404 // In case hits=0 (default), the registered hits of the current device are used.
405 // Signals which were declared as "Dead" will be rejected.
406 // The gain etc... corrected signals will be used in the process.
411 if (!hits) hits=fHits;
415 Int_t nhits=hits->GetEntries();
417 Int_t idx=0; // The signal slotindex to perform the sorting on
420 for (Int_t i=0; i<nhits; i++)
422 AliSignal* sx=(AliSignal*)hits->At(i);
426 // Obtain the slotindex corresponding to the user selection
427 idx=sx->GetSlotIndex(name);
430 if (sx->GetDeadValue(idx)) continue; // Only take alive signals
432 sig=sx->GetSignal(idx,1);
440 if (sig<vmin) vmin=sig;
441 if (sig>vmax) vmax=sig;
445 ///////////////////////////////////////////////////////////////////////////
446 TObjArray* AliDevice::SortHits(Int_t idx,Int_t mode,TObjArray* hits)
448 // Order the references to an array of hits by looping over the input array "hits"
449 // and checking the signal value. The ordered array is returned as a TObjArray.
450 // In case hits=0 (default), the registered hits of the current device are used.
451 // Note that the original hit array is not modified.
452 // A "hit" represents an abstract object which is derived from AliSignal.
453 // The user can specify the index of the signal slot to perform the sorting on.
454 // By default the slotindex will be 1.
455 // Via the "mode" argument the user can specify ordering in decreasing
456 // order (mode=-1) or ordering in increasing order (mode=1).
457 // The default is mode=-1.
458 // Signals which were declared as "Dead" will be rejected.
459 // The gain etc... corrected signals will be used in the ordering process.
467 if (!hits) hits=fHits;
469 if (idx<=0 || abs(mode)!=1 || !hits) return fOrdered;
471 Int_t nhits=hits->GetEntries();
478 fOrdered=new TObjArray(nhits);
482 for (Int_t i=0; i<nhits; i++) // Loop over all hits of the array
484 AliSignal* s=(AliSignal*)hits->At(i);
488 if (idx > s->GetNvalues()) continue; // User specified slotindex out of range for this signal
489 if (s->GetDeadValue(idx)) continue; // Only take alive signals
491 if (nord == 0) // store the first hit with a signal at the first ordered position
494 fOrdered->AddAt(s,nord-1);
498 for (Int_t j=0; j<=nord; j++) // put hit in the right ordered position
500 if (j == nord) // module has smallest (mode=-1) or largest (mode=1) signal seen so far
503 fOrdered->AddAt(s,j); // add hit at the end
504 break; // go for next hit
507 if (mode==-1 && s->GetSignal(idx,1) < ((AliSignal*)fOrdered->At(j))->GetSignal(idx,1)) continue;
508 if (mode==1 && s->GetSignal(idx,1) > ((AliSignal*)fOrdered->At(j))->GetSignal(idx,1)) continue;
511 for (Int_t k=nord-1; k>j; k--) // create empty position
513 fOrdered->AddAt(fOrdered->At(k-1),k);
515 fOrdered->AddAt(s,j); // put hit at empty position
516 break; // go for next matrix module
521 ///////////////////////////////////////////////////////////////////////////
522 TObjArray* AliDevice::SortHits(TString name,Int_t mode,TObjArray* hits)
524 // Order the references to an array of hits by looping over the input array "hits"
525 // and checking the signal value. The ordered array is returned as a TObjArray.
526 // In case hits=0 (default), the registered hits of the current device are used.
527 // Note that the input array is not modified.
528 // A "hit" represents an abstract object which is derived from AliSignal.
529 // The user can specify the name of the signal slot to perform the sorting on.
530 // In case no matching slotname is found, the signal will be skipped.
531 // Via the "mode" argument the user can specify ordering in decreasing
532 // order (mode=-1) or ordering in increasing order (mode=1).
533 // The default is mode=-1.
534 // Signals which were declared as "Dead" will be rejected.
535 // The gain etc... corrected signals will be used in the ordering process.
543 if (!hits) hits=fHits;
545 if (abs(mode)!=1 || !hits) return fOrdered;
547 Int_t nhits=hits->GetEntries();
554 fOrdered=new TObjArray(nhits);
557 Int_t idx=0; // The signal slotindex to perform the sorting on
560 for (Int_t i=0; i<nhits; i++) // loop over all hits of the array
562 AliSignal* s=(AliSignal*)hits->At(i);
566 // Obtain the slotindex corresponding to the user selection
567 idx=s->GetSlotIndex(name);
570 if (s->GetDeadValue(idx)) continue; // only take alive signals
572 if (nord == 0) // store the first hit with a signal at the first ordered position
575 fOrdered->AddAt(s,nord-1);
579 for (Int_t j=0; j<=nord; j++) // put hit in the right ordered position
581 if (j == nord) // module has smallest (mode=-1) or largest (mode=1) signal seen so far
584 fOrdered->AddAt(s,j); // add hit at the end
585 break; // go for next hit
588 if (mode==-1 && s->GetSignal(idx,1) < ((AliSignal*)fOrdered->At(j))->GetSignal(idx,1)) continue;
589 if (mode==1 && s->GetSignal(idx,1) > ((AliSignal*)fOrdered->At(j))->GetSignal(idx,1)) continue;
592 for (Int_t k=nord-1; k>j; k--) // create empty position
594 fOrdered->AddAt(fOrdered->At(k-1),k);
596 fOrdered->AddAt(s,j); // put hit at empty position
597 break; // go for next matrix module
602 ///////////////////////////////////////////////////////////////////////////
603 void AliDevice::DisplayHits(Int_t idx,Float_t scale,TObjArray* hits,Int_t dp,Int_t mstyle,Int_t mcol)
605 // 3D color display of an array hits.
606 // The user can specify the index (default=1) of the signal slot to perform the display for.
607 // The marker size will indicate the absolute value of the signal (specified by the slotindex)
608 // as a percentage of the input argument "scale".
609 // In case scale<0 the maximum absolute signal value encountered in the hit array will be used
610 // to define the 100% scale. The default is scale=-1.
611 // In case hits=0 (default), the registered hits of the current device are used.
612 // Note that the input array is not modified.
613 // In case dp=1 the device position will be used, otherwise the hit position will
614 // be used in the display. The default is dp=0.
615 // Via the "mstyle" and "mcol" arguments the user can specify the marker style
616 // and color (see TPolyMarker3D) respectively.
617 // The defaults are mstyle="large scalable dot" and mcol=blue.
618 // Signals which were declared as "Dead" will not be displayed.
619 // The gain etc... corrected signals will be used to determine the marker size.
623 // Before any display activity, a TCanvas and a TView have to be initiated
624 // first by the user like for instance
626 // TCanvas* c1=new TCanvas("c1","c1");
627 // TView* view=new TView(1);
628 // view->SetRange(-1000,-1000,-1000,1000,1000,1000);
631 Int_t thisdev=0; // Indicate whether this is the owning device or not
638 if (idx<=0 || !hits) return;
640 Int_t nhits=hits->GetEntries();
643 Float_t sigmax=fabs(scale);
647 GetExtremes(vmin,vmax,idx,hits);
649 if (fabs(vmin)>sigmax) sigmax=fabs(vmin);
652 if (sigmax <=0) return;
659 fMarkers=new TObjArray(nhits);
660 fMarkers->SetOwner();
663 GetPosition(pos,"car");
666 for (Int_t ih=0; ih<nhits; ih++)
668 AliSignal* sx=(AliSignal*)hits->At(ih);
672 sx->GetPosition(pos,"car");
678 AliDevice* dev=sx->GetDevice();
679 if (dev) dev->GetPosition(pos,"car");
682 sig=sx->GetSignal(idx,1);
683 TPolyMarker3D* m=new TPolyMarker3D();
684 m->SetMarkerStyle(mstyle);
685 m->SetMarkerColor(mcol);
686 m->SetMarkerSize(100.*fabs(sig)/sigmax);
687 m->SetPoint(0,pos[0],pos[1],pos[2]);
692 ///////////////////////////////////////////////////////////////////////////
693 void AliDevice::DisplayHits(TString name,Float_t scale,TObjArray* hits,Int_t dp,Int_t mstyle,Int_t mcol)
695 // 3D color display of an array hits.
696 // The user can specify the name of the signal slot to perform the display for.
697 // The marker size will indicate the absolute value of the signal (specified by the slotname)
698 // as a percentage of the input argument "scale".
699 // In case scale<0 the maximum absolute signal value encountered in the hit array will be used
700 // to define the 100% scale. The default is scale=-1.
701 // In case hits=0 (default), the registered hits of the current device are used.
702 // Note that the input array is not modified.
703 // In case dp=1 the device position will be used, otherwise the hit position will
704 // be used in the display. The default is dp=0.
705 // The marker size will indicate the percentage of the maximum encountered value
706 // of the absolute value of the name-specified input signal slots.
707 // Via the "mstyle" and "mcol" arguments the user can specify the marker style
708 // and color (see TPolyMarker3D) respectively.
709 // The defaults are mstyle="large scalable dot" and mcol=blue.
710 // Signals which were declared as "Dead" will not be displayed.
711 // The gain etc... corrected signals will be used to determine the marker size.
715 // Before any display activity, a TCanvas and a TView have to be initiated
716 // first by the user like for instance
718 // TCanvas* c1=new TCanvas("c1","c1");
719 // TView* view=new TView(1);
720 // view->SetRange(-1000,-1000,-1000,1000,1000,1000);
723 Int_t thisdev=0; // Indicate whether this is the owning device or not
732 Int_t nhits=hits->GetEntries();
736 Float_t sigmax=fabs(scale);
740 GetExtremes(vmin,vmax,name,hits);
742 if (fabs(vmin)>sigmax) sigmax=fabs(vmin);
745 if (sigmax <=0) return;
752 fMarkers=new TObjArray(nhits);
753 fMarkers->SetOwner();
756 GetPosition(pos,"car");
758 Int_t idx=0; // The slot index corresponding to the user specified name
760 for (Int_t ih=0; ih<nhits; ih++)
762 AliSignal* sx=(AliSignal*)hits->At(ih);
764 idx=sx->GetSlotIndex(name);
768 sx->GetPosition(pos,"car");
774 AliDevice* dev=sx->GetDevice();
775 if (dev) dev->GetPosition(pos,"car");
778 sig=sx->GetSignal(idx,1);
779 TPolyMarker3D* m=new TPolyMarker3D();
780 m->SetMarkerStyle(mstyle);
781 m->SetMarkerColor(mcol);
782 m->SetMarkerSize(100.*fabs(sig)/sigmax);
783 m->SetPoint(0,pos[0],pos[1],pos[2]);
788 ///////////////////////////////////////////////////////////////////////////
789 TObject* AliDevice::Clone(const char* name) const
791 // Make a deep copy of the current object and provide the pointer to the copy.
792 // This memberfunction enables automatic creation of new objects of the
793 // correct type depending on the object type, a feature which may be very useful
794 // for containers like AliEvent when adding objects in case the
795 // container owns the objects. This feature allows e.g. AliEvent
796 // to store either AliDevice objects or objects derived from AliDevice
797 // via tha AddDevice memberfunction, provided these derived classes also have
798 // a proper Clone memberfunction.
800 AliDevice* dev=new AliDevice(*this);
803 if (strlen(name)) dev->SetName(name);
807 ///////////////////////////////////////////////////////////////////////////