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 ///////////////////////////////////////////////////////////////////////////
19 // Class AliCalorimeter
20 // Description of a modular calorimeter system.
21 // A matrix geometry is used in which a module is identified by (row,col).
22 // Note : First module is identified as (1,1).
24 // This is the way to define and enter signals into a calorimeter :
26 // AliCalorimeter cal(10,15); // Calorimeter of 10x15 modules
27 // // All module signals set to 0.
28 // cal.AddSignal(5,7,85.4);
29 // cal.AddSignal(5,7,25.9);
30 // cal.AddSignal(3,5,1000);
31 // cal.SetSignal(5,7,10.3);
32 // cal.Reset(3,5); // Reset module (3,5) as being 'not fired'
33 // // All module data are re-initialised.
34 // cal.SetEdgeOn(1,1); // Declare module (1,1) as an 'edge module'
36 // cal.SetGain(2,8,3.2);
38 // Float_t vec[3]={6,1,20};
39 // cal.SetPosition(2,8,vec,"car");
41 // Float_t loc[3]={-1,12,3};
42 // cal.AddVetoSignal(loc,"car"); // Associate (extrapolated) position as a veto
44 // cal.Group(2); // Group 'fired' modules into clusters
45 // // Perform grouping over 2 rings around the center
46 // cal.Reset(); // Reset the complete calorimeter
47 // // Normally to prepare for the next event data
48 // // Note : Module gain, edge and dead flags remain
50 //--- Author: Nick van Eijndhoven 13-jun-1997 UU-SAP Utrecht
51 //- Modified: NvE $Date$ UU-SAP Utrecht
52 ///////////////////////////////////////////////////////////////////////////
54 #include "AliCalorimeter.h"
56 ClassImp(AliCalorimeter) // Class implementation to enable ROOT I/O
58 AliCalorimeter::AliCalorimeter()
60 // Default constructor, all parameters set to 0
77 ///////////////////////////////////////////////////////////////////////////
78 AliCalorimeter::~AliCalorimeter()
80 // Destructor to delete memory allocated to the various arrays and matrix
109 if (fMatrix || fPositions)
111 for (Int_t i=0; i<fNrows; i++)
113 for (Int_t j=0; j<fNcolumns; j++)
116 if (fPositions[i][j]) delete fPositions[i][j];
126 delete [] fPositions;
141 ///////////////////////////////////////////////////////////////////////////
142 AliCalorimeter::AliCalorimeter(Int_t nrow,Int_t ncol)
144 // Create a calorimeter module matrix
150 fAttributes=new TMatrix(nrow,ncol);
151 fGains=new TMatrix(nrow,ncol);
152 fMatrix=new AliCalmodule**[nrow];
153 fPositions=new AliPosition**[nrow];
154 for (Int_t row=0; row<nrow; row++)
156 fMatrix[row]=new AliCalmodule*[ncol];
157 fPositions[row]=new AliPosition*[ncol];
158 // Initialise the various matrices
159 for (Int_t col=0; col<ncol; col++)
162 fPositions[row][col]=0;
163 (*fGains)(row,col)=1;
164 (*fAttributes)(row,col)=0;
168 // Mark the edge modules
169 for (Int_t j=0; j<ncol; j++)
171 (*fAttributes)(0,j)=10;
172 (*fAttributes)(nrow-1,j)=10;
174 for (Int_t i=0; i<nrow; i++)
176 (*fAttributes)(i,0)=10;
177 (*fAttributes)(i,ncol-1)=10;
180 fModules=new TObjArray(); // Default size, expanded automatically
190 ///////////////////////////////////////////////////////////////////////////
191 Int_t AliCalorimeter::GetNrows()
193 // Provide the number of rows for the calorimeter module matrix
196 ///////////////////////////////////////////////////////////////////////////
197 Int_t AliCalorimeter::GetNcolumns()
199 // Provide the number of columns for the calorimeter module matrix
202 ///////////////////////////////////////////////////////////////////////////
203 void AliCalorimeter::SetSignal(Int_t row,Int_t col,Float_t sig)
205 // Set the signal for a certain calorimeter module
207 if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
209 if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
211 AliCalmodule* m=fMatrix[row-1][col-1];
212 if (!m) // only count new modules
216 AliPosition* r=fPositions[row-1][col-1];
217 if (r) m->SetPosition(*r);
219 fMatrix[row-1][col-1]=m;
221 m->SetSignal(row,col,sig);
225 cout << " *AliCalorimeter::SetSignal* row,col : " << row << "," << col
226 << " out of range." << endl;
227 cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
230 ///////////////////////////////////////////////////////////////////////////
231 void AliCalorimeter::AddSignal(Int_t row, Int_t col, Float_t sig)
233 // Add the signal to a certain calorimeter module
235 if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
237 if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
239 AliCalmodule* m=fMatrix[row-1][col-1];
240 if (!m) // only count new modules
242 SetSignal(row,col,sig);
246 m->AddSignal(row,col,sig);
251 cout << " *AliCalorimeter::AddSignal* row,col : " << row << "," << col
252 << " out of range." << endl;
253 cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
256 ///////////////////////////////////////////////////////////////////////////
257 void AliCalorimeter::AddSignal(AliCalmodule* mod)
259 // Add the signal of module mod to the current calorimeter data.
260 // This enables mixing of calorimeter data of various events.
262 if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
264 Int_t row=mod->GetRow();
265 Int_t col=mod->GetColumn();
266 Float_t sig=mod->GetSignal();
267 AliPosition r=mod->GetPosition();
269 if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
271 AliCalmodule* m=fMatrix[row-1][col-1];
272 if (!m) // No module existed yet at this position
277 fMatrix[row-1][col-1]=m;
280 m->AddSignal(row,col,sig);
281 if (!fPositions[row-1][col-1]) fPositions[row-1][col-1]=new AliPosition(r);
285 cout << " *AliCalorimeter::AddSignal* row,col : " << row << "," << col
286 << " out of range." << endl;
287 cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
290 ///////////////////////////////////////////////////////////////////////////
291 void AliCalorimeter::Reset(Int_t row,Int_t col)
293 // Reset the signal for a certain calorimeter module
295 if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
297 if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
299 AliCalmodule* m=fMatrix[row-1][col-1];
304 fModules->Compress();
306 fMatrix[row-1][col-1]=0;
311 cout << " *AliCalorimeter::Reset* row,col : " << row << "," << col
312 << " out of range." << endl;
313 cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
316 ///////////////////////////////////////////////////////////////////////////
317 void AliCalorimeter::Reset()
319 // Reset the signals for the complete calorimeter
320 // Normally this is done to prepare for the data of the next event
321 // Note : Module gains, edge and dead flags remain unchanged
323 if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
326 if (fModules) fModules->Delete();
327 for (Int_t i=0; i<fNrows; i++)
329 for (Int_t j=0; j<fNcolumns; j++)
351 ///////////////////////////////////////////////////////////////////////////
352 Float_t AliCalorimeter::GetSignal(Int_t row,Int_t col)
354 // Provide the signal of a certain calorimeter module.
355 // In case the module was marked dead, 0 is returned.
357 if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
359 if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
362 AliCalmodule* m=fMatrix[row-1][col-1];
365 Int_t dead=m->GetDeadValue();
366 if (!dead) signal=m->GetSignal();
372 cout << " *AliCalorimeter::GetSignal* row,col : " << row << "," << col
373 << " out of range." << endl;
374 cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
378 ///////////////////////////////////////////////////////////////////////////
379 void AliCalorimeter::SetEdgeOn(Int_t row,Int_t col)
381 // Indicate a certain calorimeter module as 'edge module'
383 if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
385 if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
387 Float_t word=(*fAttributes)(row-1,col-1);
388 Int_t iword=int(word+0.1);
391 (*fAttributes)(row-1,col-1)=float(dead+10*edge);
395 cout << " *AliCalorimeter::SetEdgeOn* row,col : " << row << "," << col
396 << " out of range." << endl;
397 cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
400 ///////////////////////////////////////////////////////////////////////////
401 void AliCalorimeter::SetEdgeOff(Int_t row,Int_t col)
403 // Indicate a certain calorimeter module as 'non-edge module'
405 if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
407 if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
409 Float_t word=(*fAttributes)(row-1,col-1);
410 Int_t iword=int(word+0.1);
413 (*fAttributes)(row-1,col-1)=float(dead+10*edge);
417 cout << " *AliCalorimeter::SetEdgeOff* row,col : " << row << "," << col
418 << " out of range." << endl;
419 cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
422 ///////////////////////////////////////////////////////////////////////////
423 void AliCalorimeter::SetDead(Int_t row,Int_t col)
425 // Indicate a certain calorimeter module as 'dead module'
427 if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
429 if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
431 Float_t word=(*fAttributes)(row-1,col-1);
432 Int_t iword=int(word+0.1);
435 (*fAttributes)(row-1,col-1)=float(dead+10*edge);
436 if (fMatrix[row-1][col-1]) (fMatrix[row-1][col-1])->SetDead();
438 // Increase the 'edge value' of surrounding modules
444 if (rlow < 1) rlow=row;
445 if (rup > fNrows) rup=fNrows;
446 if (clow < 1) clow=col;
447 if (cup > fNcolumns) cup=fNcolumns;
449 for (Int_t i=rlow; i<=rup; i++)
451 for (Int_t j=clow; j<=cup; j++)
453 if (i!=row || j!=col) // No increase of edge value for the dead module itself
455 word=(*fAttributes)(i-1,j-1);
460 (*fAttributes)(i-1,j-1)=float(dead+10*edge);
467 cout << " *AliCalorimeter::SetDead* row,col : " << row << "," << col
468 << " out of range." << endl;
469 cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
472 ///////////////////////////////////////////////////////////////////////////
473 void AliCalorimeter::SetAlive(Int_t row,Int_t col)
475 // Indicate a certain calorimeter module as 'active module'
477 if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
479 if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
481 Float_t word=(*fAttributes)(row-1,col-1);
482 Int_t iword=int(word+0.1);
485 (*fAttributes)(row-1,col-1)=float(dead+10*edge);
486 if (fMatrix[row-1][col-1]) (fMatrix[row-1][col-1])->SetAlive();
488 // Decrease the 'edge value' of surrounding modules
494 if (rlow < 1) rlow=row;
495 if (rup > fNrows) rup=fNrows;
496 if (clow < 1) clow=col;
497 if (cup > fNcolumns) cup=fNcolumns;
499 for (Int_t i=rlow; i<=rup; i++)
501 for (Int_t j=clow; j<=cup; j++)
503 if (i!=row || j!=col) // No decrease of edge value for the dead module itself
505 word=(*fAttributes)(i-1,j-1);
510 (*fAttributes)(i-1,j-1)=float(dead+10*edge);
517 cout << " *AliCalorimeter::SetAlive* row,col : " << row << "," << col
518 << " out of range." << endl;
519 cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
522 ///////////////////////////////////////////////////////////////////////////
523 void AliCalorimeter::SetGain(Int_t row,Int_t col,Float_t gain)
525 // Set the gain value for a certain calorimeter module
527 if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
529 if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
531 (*fGains)(row-1,col-1)=gain;
532 if (fMatrix[row-1][col-1]) (fMatrix[row-1][col-1])->SetGain(gain);
536 cout << " *AliCalorimeter::SetGain* row,col : " << row << "," << col
537 << " out of range." << endl;
538 cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
541 ///////////////////////////////////////////////////////////////////////////
542 void AliCalorimeter::SetPosition(Int_t row,Int_t col,Float_t* vec,TString f)
544 // Set the position in user coordinates for a certain calorimeter module
546 if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
548 if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
550 if (!fPositions[row-1][col-1]) fPositions[row-1][col-1]=new AliPosition;
551 (fPositions[row-1][col-1])->SetVector(vec,f);
552 if (fMatrix[row-1][col-1]) (fMatrix[row-1][col-1])->SetPosition(vec,f);
556 cout << " *AliCalorimeter::SetPosition* row,col : " << row << "," << col
557 << " out of range." << endl;
558 cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
561 ///////////////////////////////////////////////////////////////////////////
562 void AliCalorimeter::SetPosition(Int_t row,Int_t col,Ali3Vector& r)
564 // Set the position for a certain calorimeter module
566 if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
568 if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
570 if (!fPositions[row-1][col-1]) fPositions[row-1][col-1]=new AliPosition;
571 (fPositions[row-1][col-1])->SetPosition(r);
572 if (fMatrix[row-1][col-1]) (fMatrix[row-1][col-1])->SetPosition(r);
576 cout << " *AliCalorimeter::SetPosition* row,col : " << row << "," << col
577 << " out of range." << endl;
578 cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
581 ///////////////////////////////////////////////////////////////////////////
582 Int_t AliCalorimeter::GetEdgeValue(Int_t row,Int_t col)
584 // Provide the value of the edge flag of a certain module
586 if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
588 Float_t word=(*fAttributes)(row-1,col-1);
589 Int_t iword=int(word+0.1);
595 cout << " *AliCalorimeter::GetEdgeValue* row,col : " << row << "," << col
596 << " out of range." << endl;
597 cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
601 ///////////////////////////////////////////////////////////////////////////
602 Int_t AliCalorimeter::GetDeadValue(Int_t row,Int_t col)
604 // Provide the value of the dead flag of a certain module
606 if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
608 Float_t word=(*fAttributes)(row-1,col-1);
609 Int_t iword=int(word+0.1);
615 cout << " *AliCalorimeter::GetDeadValue* row,col : " << row << "," << col
616 << " out of range." << endl;
617 cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
621 ///////////////////////////////////////////////////////////////////////////
622 Float_t AliCalorimeter::GetGain(Int_t row,Int_t col)
624 // Provide the gain value of a certain module
626 if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
628 return (*fGains)(row-1,col-1);
632 cout << " *AliCalorimeter::GetGain* row,col : " << row << "," << col
633 << " out of range." << endl;
634 cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
638 ///////////////////////////////////////////////////////////////////////////
639 void AliCalorimeter::GetPosition(Int_t row,Int_t col,Float_t* vec,TString f)
641 // Return the position in user coordinates for a certain calorimeter module
643 if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
645 if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
647 // if (fMatrix[row-1][col-1]) (fMatrix[row-1][col-1])->GetPosition(vec,f);
648 if (fPositions[row-1][col-1]) (fPositions[row-1][col-1])->GetVector(vec,f);
652 cout << " *AliCalorimeter::GetPosition* row,col : " << row << "," << col
653 << " out of range." << endl;
654 cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
657 ///////////////////////////////////////////////////////////////////////////
658 AliPosition* AliCalorimeter::GetPosition(Int_t row,Int_t col)
660 // Access to the position of a certain calorimeter module
662 if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
664 if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
666 return fPositions[row-1][col-1];
670 cout << " *AliCalorimeter::GetPosition* row,col : " << row << "," << col
671 << " out of range." << endl;
672 cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
676 ///////////////////////////////////////////////////////////////////////////
677 Float_t AliCalorimeter::GetClusteredSignal(Int_t row,Int_t col)
679 // Provide the module signal after clustering
681 if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
683 if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
685 if (fMatrix[row-1][col-1])
687 return (fMatrix[row-1][col-1])->GetClusteredSignal();
696 cout << " *AliCalorimeter::GetClusteredSignal* row,col : " << row << "," << col
697 << " out of range." << endl;
698 cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
702 ///////////////////////////////////////////////////////////////////////////
703 Int_t AliCalorimeter::GetNsignals()
705 // Provide the number of modules that contain a signal
706 // Note : The number of modules marked 'dead' but which had a signal
710 ///////////////////////////////////////////////////////////////////////////
711 void AliCalorimeter::Group(Int_t n)
713 // Group the individual modules into clusters
714 // Module signals of n rings around the central module will be grouped
716 if (fNsignals > 0) // Directly return if no modules fired
718 if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
720 if (fNclusters > 0) Ungroup(); // Restore unclustered situation if needed
722 // Order the modules with decreasing signal
723 AliCalmodule** ordered=new AliCalmodule*[fNsignals]; // temp. array for ordered modules
727 // Clustering of modules. Start with the highest signal.
734 fClusters=new TObjArray();
740 for (Int_t i=0; i<nord; i++)
742 row=ordered[i]->GetRow(); // row number of cluster center
743 col=ordered[i]->GetColumn(); // column number of cluster center
744 if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
746 m=fMatrix[row-1][col-1];
749 // only use modules not yet used in a cluster
750 if (m->GetClusteredSignal() > 0.)
752 Int_t edge=GetEdgeValue(row,col);
754 if (!edge) c->Start(*m); // module to start the cluster if not on edge
755 if (c->GetNmodules() > 0) // cluster started successfully (no edge)
758 fNclusters++; // update cluster counter
759 AddRing(row,col,n); // add signals of n rings around the center
770 // Delete the temp. array
773 for (Int_t j=0; j<nord; j++)
781 ///////////////////////////////////////////////////////////////////////////
782 void AliCalorimeter::Sortm(AliCalmodule** ordered,Int_t& nord)
784 // Order the modules with decreasing signal
787 for (Int_t i=0; i<fNrows; i++) // loop over all modules of the matrix
789 for (Int_t ii=0; ii<fNcolumns; ii++)
791 if (GetSignal(i+1,ii+1) <= 0.) continue; // only take alive modules with a signal
793 if (nord == 0) // store the first module with a signal at the first ordered position
796 ordered[nord-1]=fMatrix[i][ii];
800 for (Int_t j=0; j<=nord; j++) // put module in the right ordered position
802 if (j == nord) // module has smallest signal seen so far
805 ordered[j]=fMatrix[i][ii]; // add module at the end
806 break; // go for next matrix module
809 if (GetSignal(i+1,ii+1) < ordered[j]->GetSignal()) continue;
812 for (Int_t k=nord-1; k>j; k--) // create empty position
814 ordered[k]=ordered[k-1];
816 ordered[j]=fMatrix[i][ii]; // put module at empty position
817 break; // go for next matrix module
822 ///////////////////////////////////////////////////////////////////////////
823 void AliCalorimeter::AddRing(Int_t row, Int_t col, Int_t n)
825 // Add module signals of 1 ring around (row,col) to current cluster
826 // n denotes the maximum number of rings around cluster center
827 // Note : This function is used recursively
829 if (n >= 1) // Check if any rings left for recursive calls
831 Float_t signal=GetSignal(row,col); // signal of (row,col) module
833 Int_t lrow=row-1; if (lrow < 1) lrow=1; // row lowerbound for ring
834 Int_t urow=row+1; if (urow > fNrows) urow=fNrows; // row upperbound for ring
835 Int_t lcol=col-1; if (lcol < 1) lcol=1; // col lowerbound for ring
836 Int_t ucol=col+1; if (ucol > fNcolumns) ucol=fNcolumns; // row upperbound for ring
838 for (Int_t i=lrow; i<=urow; i++)
840 for (Int_t j=lcol; j<=ucol; j++)
842 // add module(i,j) to cluster if the signal <= signal(row,col)
843 if (GetSignal(i,j) <= signal)
845 AliCalmodule* m=fMatrix[i-1][j-1];
846 if (m) ((AliCalcluster*)fClusters->At(fNclusters-1))->Add(*m);
848 AddRing(i,j,n-1); // Go for ring of modules around this (i,j) one
853 ///////////////////////////////////////////////////////////////////////////
854 Int_t AliCalorimeter::GetNclusters()
856 // Provide the number of clusters
859 ///////////////////////////////////////////////////////////////////////////
860 AliCalcluster* AliCalorimeter::GetCluster(Int_t j)
862 // Provide cluster number j
863 // Note : j=1 denotes the first cluster
864 if ((j >= 1) && (j <= fNclusters))
866 return (AliCalcluster*)fClusters->At(j-1);
870 cout << " *AliCalorimeter::GetCluster* cluster number : " << j
871 << " out of range." << endl;
872 cout << " -- Cluster number 1 (if any) returned " << endl;
873 return (AliCalcluster*)fClusters->At(0);
876 ///////////////////////////////////////////////////////////////////////////
877 AliCalmodule* AliCalorimeter::GetModule(Int_t j)
879 // Provide 'fired' module number j
880 // Note : j=1 denotes the first 'fired' module
881 if ((j >= 1) && (j <= fNsignals))
883 return (AliCalmodule*)fModules->At(j-1);
887 cout << " *AliCalorimeter::GetModule* module number : " << j
888 << " out of range." << endl;
889 cout << " -- Fired module number 1 (if any) returned " << endl;
890 return (AliCalmodule*)fModules->At(0);
893 ///////////////////////////////////////////////////////////////////////////
894 AliCalmodule* AliCalorimeter::GetModule(Int_t row,Int_t col)
896 // Provide access to module (row,col).
897 // Note : first module is at (1,1).
899 if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
901 if (row>=1 && row<=fNrows && col>=1 && col<=fNcolumns)
903 return fMatrix[row-1][col-1];
907 cout << " *AliCalorimeter::GetModule* row,col : " << row << ", " << col
908 << " out of range." << endl;
912 ///////////////////////////////////////////////////////////////////////////
913 TH2F* AliCalorimeter::DrawModules()
915 // Provide a lego plot of the module signals
923 fHmodules=new TH2F("fHmodules","Module signals",
924 fNcolumns,0.5,float(fNcolumns)+0.5,fNrows,0.5,float(fNrows)+0.5);
926 fHmodules->SetDirectory(0); // Suppress global character of histo pointer
930 Float_t row,col,signal;
932 for (Int_t i=0; i<fNsignals; i++)
934 m=(AliCalmodule*)fModules->At(i);
937 row=float(m->GetRow());
938 col=float(m->GetColumn());
939 dead=m->GetDeadValue();
941 if (!dead) signal=m->GetSignal();
942 if (signal>0.) fHmodules->Fill(col,row,signal);
946 fHmodules->Draw("lego");
949 ///////////////////////////////////////////////////////////////////////////
950 TH2F* AliCalorimeter::DrawClusters()
952 // Provide a lego plot of the cluster signals
960 fHclusters=new TH2F("fHclusters","Cluster signals",
961 fNcolumns,0.5,float(fNcolumns)+0.5,fNrows,0.5,float(fNrows)+0.5);
963 fHclusters->SetDirectory(0); // Suppress global character of histo pointer
967 Float_t row,col,signal;
968 for (Int_t i=0; i<fNclusters; i++)
970 c=(AliCalcluster*)fClusters->At(i);
973 row=float(c->GetRow());
974 col=float(c->GetColumn());
975 signal=c->GetSignal();
976 if (signal>0.) fHclusters->Fill(col,row,signal);
980 fHclusters->Draw("lego");
983 ///////////////////////////////////////////////////////////////////////////
984 void AliCalorimeter::LoadMatrix()
986 // Load the Calorimeter module matrix data back from the TObjArray
988 // Initialise the module matrix
991 fMatrix=new AliCalmodule**[fNrows];
992 for (Int_t i=0; i<fNrows; i++)
994 fMatrix[i]=new AliCalmodule*[fNcolumns];
998 // Initialise the position matrix
1001 fPositions=new AliPosition**[fNrows];
1002 for (Int_t j=0; j<fNrows; j++)
1004 fPositions[j]=new AliPosition*[fNcolumns];
1008 for (Int_t jrow=0; jrow<fNrows; jrow++)
1010 for (Int_t jcol=0; jcol<fNcolumns; jcol++)
1012 fMatrix[jrow][jcol]=0;
1013 fPositions[jrow][jcol]=0;
1017 // Copy the module pointers back into the matrix
1022 if (fModules) nsig=fModules->GetEntries();
1023 for (Int_t j=0; j<nsig; j++)
1025 m=(AliCalmodule*)fModules->At(j);
1030 AliPosition r=m->GetPosition();
1031 fMatrix[row-1][col-1]=m;
1032 fPositions[row-1][col-1]=new AliPosition(r);
1036 ///////////////////////////////////////////////////////////////////////////
1037 void AliCalorimeter::Ungroup()
1039 // Set the module signals back to the non-clustered situation
1041 if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
1044 if (fModules) nsig=fModules->GetEntries();
1048 for (Int_t j=0; j<nsig; j++)
1050 m=(AliCalmodule*)fModules->At(j);
1053 signal=m->GetSignal();
1054 m->SetClusteredSignal(signal);
1058 ///////////////////////////////////////////////////////////////////////////
1059 void AliCalorimeter::AddVetoSignal(Float_t* r,TString f,Float_t s)
1061 // Associate an (extrapolated) AliSignal at location r as veto to the cal.
1062 // Note : The default signal value (s) is 0
1066 fVetos=new TObjArray();
1069 fVetos->Add(new AliSignal);
1072 ((AliSignal*)fVetos->At(fNvetos-1))->SetPosition(r,f);
1073 ((AliSignal*)fVetos->At(fNvetos-1))->SetSignal(s);
1075 ///////////////////////////////////////////////////////////////////////////
1076 Int_t AliCalorimeter::GetNvetos()
1078 // Provide the number of veto signals associated to the calorimeter
1081 ///////////////////////////////////////////////////////////////////////////
1082 AliSignal* AliCalorimeter::GetVetoSignal(Int_t i)
1084 // Provide access to the i-th veto signal of this calorimeter
1085 // Note : The first hit corresponds to i=1
1087 if (i>0 && i<=fNvetos)
1089 return (AliSignal*)fVetos->At(i-1);
1093 cout << " *AliCalorimeter::GetVetoSignal* Signal number " << i
1094 << " out of range." << endl;
1095 cout << " --- First signal (if any) returned." << endl;
1096 return (AliSignal*)fVetos->At(0);
1099 ///////////////////////////////////////////////////////////////////////////
1100 void AliCalorimeter::SetName(TString name)
1104 ///////////////////////////////////////////////////////////////////////////
1105 TString AliCalorimeter::GetName()
1109 ///////////////////////////////////////////////////////////////////////////