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.
30 // // Set user defined status word to indicate e.g. readout electronics version
31 // m.SetStatus(100201);
33 // m.SetName("OM123");
35 // Float_t pos[3]={1,2,3};
36 // m.SetPosition(pos,"car");
41 // s.SetName("OM123 Hit 1");
42 // s.SetSlotName("ADC");
44 // s.SetSlotName("LE",2);
45 // s.SetSignal(-100,2);
46 // s.SetSlotName("TOT",3);
47 // s.SetSignal(-1000,3);
51 // s.SetName("OM123 Hit 2");
52 // s.SetSlotName("ADC");
54 // s.SetSlotName("LE",2);
55 // s.SetSignal(-101,2);
56 // s.SetSlotName("TOT",3);
57 // s.SetSignal(1001,3);
61 // s.SetName("OM123 Hit 3");
62 // s.SetSlotName("ADC");
64 // s.SetSlotName("LE",2);
65 // s.SetSignal(-102,2);
66 // s.SetSlotName("TOT",3);
67 // s.SetSignal(-1002,3);
70 // TObjArray* ordered=m.SortHits("TOT");
71 // nhits=ordered->GetEntries();
72 // for (Int_t i=0; i<nhits; i++)
74 // AliSignal* sx=(AliSignal*)ordered->At(i);
75 // if (sx) sx->Data();
78 //--- Author: Nick van Eijndhoven 23-jun-2004 Utrecht University
79 //- Modified: NvE $Date$ Utrecht University
80 ///////////////////////////////////////////////////////////////////////////
83 #include "AliDevice.h"
84 #include "Riostream.h"
86 ClassImp(AliDevice) // Class implementation to enable ROOT I/O
88 AliDevice::AliDevice() : AliSignal()
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.
101 ///////////////////////////////////////////////////////////////////////////
102 AliDevice::~AliDevice()
104 // Default destructor.
106 // Remove backward links to this device from the hits
107 // which were not owned by it.
110 for (Int_t ih=1; ih<=GetNhits(); ih++)
112 AliSignal* sx=GetHit(ih);
113 if (sx) sx->ResetLinks(this);
135 ///////////////////////////////////////////////////////////////////////////
136 AliDevice::AliDevice(const AliDevice& dev) : AliSignal(dev)
144 fStatus=dev.GetStatus();
145 fHitCopy=dev.GetHitCopy();
147 Int_t nhits=dev.GetNhits();
150 fHits=new TObjArray(nhits);
151 if (fHitCopy) fHits->SetOwner();
152 for (Int_t ih=1; ih<=nhits; ih++)
154 AliSignal* sx=dev.GetHit(ih);
157 fHits->Add(sx->Clone());
158 AliSignal* s=(AliSignal*)fHits->Last();
159 s->ResetLinks((AliDevice*)&dev);
170 ///////////////////////////////////////////////////////////////////////////
171 void AliDevice::Reset(Int_t mode)
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.
178 AliSignal::Reset(mode);
180 ///////////////////////////////////////////////////////////////////////////
181 void AliDevice::SetHitCopy(Int_t j)
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.
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().
199 cout << "*AliDevice::SetHitCopy* Invalid argument : " << j << endl;
204 cout << "*AliDevice::SetHitCopy* Storage already contained hits."
205 << " ==> HitCopy mode not changed." << endl;
208 ///////////////////////////////////////////////////////////////////////////
209 Int_t AliDevice::GetHitCopy() const
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.
216 ///////////////////////////////////////////////////////////////////////////
217 void AliDevice::SetStatus(Int_t word)
219 // Set a user defined status word for this device.
222 ///////////////////////////////////////////////////////////////////////////
223 Int_t AliDevice::GetStatus() const
225 // Provide the user defined status word for this device.
228 ///////////////////////////////////////////////////////////////////////////
229 void AliDevice::AddHit(AliSignal& s)
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.
241 fHits=new TObjArray(1);
242 if (fHitCopy) fHits->SetOwner();
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++)
249 if (&s==fHits->At(i)) return;
252 // Check for existing (backward) link to this device.
253 Int_t nlinks=s.GetNlinks(this);
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);
267 // Set (backward) link to the this device
268 if (!nlinks) s.AddLink(this);
271 ///////////////////////////////////////////////////////////////////////////
272 void AliDevice::RemoveHit(AliSignal& s)
274 // Remove AliSignal object registered as a hit from this device.
277 AliSignal* test=(AliSignal*)fHits->Remove(&s);
281 if (fHitCopy) delete test;
286 AliSignal* test=(AliSignal*)fOrdered->Remove(&s);
287 if (test) fOrdered->Compress();
290 ///////////////////////////////////////////////////////////////////////////
291 void AliDevice::RemoveHits()
293 // Remove all AliSignal objects registered as hits from this device.
310 ///////////////////////////////////////////////////////////////////////////
311 Int_t AliDevice::GetNhits() const
313 // Provide the number of registered hits for this device.
315 if (fHits) nhits=fHits->GetEntries();
318 ///////////////////////////////////////////////////////////////////////////
319 AliSignal* AliDevice::GetHit(Int_t j) const
321 // Provide the AliSignal object registered as hit number j.
322 // Note : j=1 denotes the first hit.
323 if (!fHits) return 0;
325 if ((j >= 1) && (j <= GetNhits()))
327 return (AliSignal*)fHits->At(j-1);
334 ///////////////////////////////////////////////////////////////////////////
335 AliSignal* AliDevice::GetHit(TString name) const
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.
340 if (!fHits) return 0;
343 Int_t nhits=GetNhits();
344 for (Int_t i=0; i<nhits; i++)
346 AliSignal* sx=(AliSignal*)fHits->At(i);
349 hitname=sx->GetName();
350 if (hitname == name) return sx;
354 return 0; // No matching name found
356 ///////////////////////////////////////////////////////////////////////////
357 AliSignal* AliDevice::GetIdHit(Int_t id) const
359 // Return the hit with unique identifier "id".
360 if (!fHits || id<0) return 0;
364 for (Int_t i=0; i<GetNhits(); i++)
366 sx=(AliSignal*)fHits->At(i);
369 sid=sx->GetUniqueID();
370 if (id==sid) return sx;
373 return 0; // No matching id found
375 ///////////////////////////////////////////////////////////////////////////
376 TObjArray* AliDevice::GetHits()
378 // Provide the references to all the registered hits.
381 ///////////////////////////////////////////////////////////////////////////
382 void AliDevice::ShowHit(Int_t j) const
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.
389 Int_t nhits=GetNhits();
390 for (Int_t ih=1; ih<=nhits; ih++)
392 AliSignal* sx=GetHit(ih);
398 AliSignal* s=GetHit(j);
402 ///////////////////////////////////////////////////////////////////////////
403 void AliDevice::Data(TString f,TString u) const
405 // Print the device and all registered hit info according to the specified
406 // coordinate frame f.
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
413 // The defaults are f="car" and u="rad".
415 AliSignal::Data(f,u);
416 Int_t nhits=GetNhits();
419 cout << " The following " << nhits << " hits are registered : " << endl;
424 cout << " No hits have been registered for this device." << endl;
427 ///////////////////////////////////////////////////////////////////////////
428 void AliDevice::GetExtremes(Float_t& vmin,Float_t& vmax,Int_t idx,TObjArray* hits,Int_t mode) const
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).
443 if (!hits) hits=fHits;
445 if (idx<=0 || !hits) return;
447 Int_t nhits=hits->GetEntries();
450 for (Int_t i=0; i<nhits; i++)
452 AliSignal* sx=(AliSignal*)hits->At(i);
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
458 sig=sx->GetSignal(idx,mode);
466 if (sig<vmin) vmin=sig;
467 if (sig>vmax) vmax=sig;
471 ///////////////////////////////////////////////////////////////////////////
472 void AliDevice::GetExtremes(Float_t& vmin,Float_t& vmax,TString name,TObjArray* hits,Int_t mode) const
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).
486 if (!hits) hits=fHits;
490 Int_t nhits=hits->GetEntries();
492 Int_t idx=0; // The signal slotindex to perform the sorting on
495 for (Int_t i=0; i<nhits; i++)
497 AliSignal* sx=(AliSignal*)hits->At(i);
501 // Obtain the slotindex corresponding to the user selection
502 idx=sx->GetSlotIndex(name);
505 if (sx->GetDeadValue(idx)) continue; // Only take alive signals
507 sig=sx->GetSignal(idx,mode);
515 if (sig<vmin) vmin=sig;
516 if (sig>vmax) vmax=sig;
520 ///////////////////////////////////////////////////////////////////////////
521 TObjArray* AliDevice::SortHits(Int_t idx,Int_t mode,TObjArray* hits,Int_t mcal)
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).
546 if (!hits) hits=fHits;
548 if (idx<=0 || abs(mode)!=1 || !hits) return fOrdered;
550 Int_t nhits=hits->GetEntries();
557 fOrdered=new TObjArray(nhits);
561 for (Int_t i=0; i<nhits; i++) // Loop over all hits of the array
563 AliSignal* s=(AliSignal*)hits->At(i);
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
570 if (nord == 0) // store the first hit with a signal at the first ordered position
573 fOrdered->AddAt(s,nord-1);
577 for (Int_t j=0; j<=nord; j++) // put hit in the right ordered position
579 if (j == nord) // module has smallest (mode=-1) or largest (mode=1) signal seen so far
582 fOrdered->AddAt(s,j); // add hit at the end
583 break; // go for next hit
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;
590 for (Int_t k=nord-1; k>j; k--) // create empty position
592 fOrdered->AddAt(fOrdered->At(k-1),k);
594 fOrdered->AddAt(s,j); // put hit at empty position
595 break; // go for next hit
600 ///////////////////////////////////////////////////////////////////////////
601 TObjArray* AliDevice::SortHits(TString name,Int_t mode,TObjArray* hits,Int_t mcal)
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).
626 if (!hits) hits=fHits;
628 if (abs(mode)!=1 || !hits) return fOrdered;
630 Int_t nhits=hits->GetEntries();
637 fOrdered=new TObjArray(nhits);
640 Int_t idx=0; // The signal slotindex to perform the sorting on
643 for (Int_t i=0; i<nhits; i++) // loop over all hits of the array
645 AliSignal* s=(AliSignal*)hits->At(i);
649 // Obtain the slotindex corresponding to the user selection
650 idx=s->GetSlotIndex(name);
653 if (s->GetDeadValue(idx)) continue; // only take alive signals
655 if (nord == 0) // store the first hit with a signal at the first ordered position
658 fOrdered->AddAt(s,nord-1);
662 for (Int_t j=0; j<=nord; j++) // put hit in the right ordered position
664 if (j == nord) // module has smallest (mode=-1) or largest (mode=1) signal seen so far
667 fOrdered->AddAt(s,j); // add hit at the end
668 break; // go for next hit
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;
675 for (Int_t k=nord-1; k>j; k--) // create empty position
677 fOrdered->AddAt(fOrdered->At(k-1),k);
679 fOrdered->AddAt(s,j); // put hit at empty position
680 break; // go for next hit
685 ///////////////////////////////////////////////////////////////////////////
686 void AliDevice::DisplayHits(Int_t idx,Float_t scale,TObjArray* hits,Int_t dp,Int_t mode,Int_t mcol)
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).
709 // Before any display activity, a TCanvas and a TView have to be initiated
710 // first by the user like for instance
712 // TCanvas* c1=new TCanvas("c1","c1");
713 // TView* view=new TView(1);
714 // view->SetRange(-1000,-1000,-1000,1000,1000,1000);
717 Int_t thisdev=0; // Indicate whether this is the owning device or not
724 if (idx<=0 || !hits) return;
726 Int_t nhits=hits->GetEntries();
729 Float_t sigmax=fabs(scale);
733 GetExtremes(vmin,vmax,idx,hits,mode);
735 if (fabs(vmin)>sigmax) sigmax=fabs(vmin);
738 if (sigmax <=0) return;
745 fMarkers=new TObjArray(nhits);
746 fMarkers->SetOwner();
749 GetPosition(pos,"car");
752 for (Int_t ih=0; ih<nhits; ih++)
754 AliSignal* sx=(AliSignal*)hits->At(ih);
758 sx->GetPosition(pos,"car");
764 AliDevice* dev=sx->GetDevice();
765 if (dev) dev->GetPosition(pos,"car");
768 sig=sx->GetSignal(idx,mode);
771 if (fabs(sig) <= 0.) continue;
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]);
782 ///////////////////////////////////////////////////////////////////////////
783 void AliDevice::DisplayHits(TString name,Float_t scale,TObjArray* hits,Int_t dp,Int_t mode,Int_t mcol)
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).
808 // Before any display activity, a TCanvas and a TView have to be initiated
809 // first by the user like for instance
811 // TCanvas* c1=new TCanvas("c1","c1");
812 // TView* view=new TView(1);
813 // view->SetRange(-1000,-1000,-1000,1000,1000,1000);
816 Int_t thisdev=0; // Indicate whether this is the owning device or not
825 Int_t nhits=hits->GetEntries();
829 Float_t sigmax=fabs(scale);
833 GetExtremes(vmin,vmax,name,hits,mode);
835 if (fabs(vmin)>sigmax) sigmax=fabs(vmin);
838 if (sigmax <=0) return;
845 fMarkers=new TObjArray(nhits);
846 fMarkers->SetOwner();
849 GetPosition(pos,"car");
851 Int_t idx=0; // The slot index corresponding to the user specified name
853 for (Int_t ih=0; ih<nhits; ih++)
855 AliSignal* sx=(AliSignal*)hits->At(ih);
857 idx=sx->GetSlotIndex(name);
861 sx->GetPosition(pos,"car");
867 AliDevice* dev=sx->GetDevice();
868 if (dev) dev->GetPosition(pos,"car");
871 sig=sx->GetSignal(idx,mode);
874 if (fabs(sig) <= 0.) continue;
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]);
885 ///////////////////////////////////////////////////////////////////////////
886 TObject* AliDevice::Clone(const char* name) const
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.
897 AliDevice* dev=new AliDevice(*this);
900 if (strlen(name)) dev->SetName(name);
904 ///////////////////////////////////////////////////////////////////////////