]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RALICE/AliCalorimeter.cxx
08-jan-2003 NvE Delete statements for fMatrix[i] and fPositions[i] added in the dtor of
[u/mrichter/AliRoot.git] / RALICE / AliCalorimeter.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 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).
23 //
24 // This is the way to define and enter signals into a calorimeter :
25 //
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'
35 //   cal.SetDead(8,3);
36 //   cal.SetGain(2,8,3.2);
37 //
38 //   Float_t vec[3]={6,1,20};
39 //   cal.SetPosition(2,8,vec,"car");
40 //
41 //   AliSignal s;
42 //   Float_t loc[3]={-1,12,3};
43 //   s.SetPosition(loc,"car");
44 //   s.SetSignal(328);
45 //   cal.AddVetoSignal(s); // Associate (extrapolated) signal as a veto
46 //
47 //   cal.Group(2);      // Group 'fired' modules into clusters
48 //                      // Perform grouping over 2 rings around the center
49 //   cal.Reset();       // Reset the complete calorimeter
50 //                      // Normally to prepare for the next event data
51 //                      // Note : Module gain, edge and dead flags remain
52 //
53 //--- Author: Nick van Eijndhoven 13-jun-1997 UU-SAP Utrecht
54 //- Modified: NvE $Date$ UU-SAP Utrecht
55 ///////////////////////////////////////////////////////////////////////////
56
57 #include "AliCalorimeter.h"
58
59 ClassImp(AliCalorimeter) // Class implementation to enable ROOT I/O
60  
61 AliCalorimeter::AliCalorimeter()
62 {
63 // Default constructor, all parameters set to 0
64  fNrows=0;
65  fNcolumns=0;
66  fNsignals=0;
67  fNclusters=0;
68  fMatrix=0;
69  fClusters=0;
70  fModules=0;
71  fHmodules=0;
72  fHclusters=0;
73  fNvetos=0;
74  fVetos=0;
75  fAttributes=0;
76  fGains=0;
77  fPositions=0;
78  fName=" ";
79 }
80 ///////////////////////////////////////////////////////////////////////////
81 AliCalorimeter::~AliCalorimeter()
82 {
83 // Destructor to delete memory allocated to the various arrays and matrices
84  if (fModules)
85  {
86   delete fModules;
87   fModules=0;
88  }
89  if (fClusters)
90  {
91   delete fClusters;
92   fClusters=0;
93  }
94  if (fVetos)
95  {
96   delete fVetos;
97   fVetos=0;
98  }
99  if (fHmodules)
100  {
101   delete fHmodules;
102   fHmodules=0;
103  }
104  if (fHclusters)
105  {
106   delete fHclusters;
107   fHclusters=0;
108  }
109
110  // Free memory allocated for (internal) module and position matrices.
111  // The modules have already been deleted via the fModules I/O array,
112  // so they shouldn't be deleted here anymore.
113  if (fMatrix || fPositions)
114  {
115   for (Int_t i=0; i<fNrows; i++)
116   {
117    for (Int_t j=0; j<fNcolumns; j++)
118    {
119     if (fPositions[i][j]) delete fPositions[i][j];
120    }
121    if (fPositions[i]) delete [] fPositions[i];
122    if (fMatrix[i]) delete [] fMatrix[i];
123   }
124   if (fMatrix)
125   {
126    delete [] fMatrix;
127    fMatrix=0;
128   }
129   if (fPositions)
130   {
131    delete [] fPositions;
132    fPositions=0;
133   }
134  }
135  if (fGains)
136  {
137   delete fGains;
138   fGains=0;
139  }
140  if (fAttributes)
141  {
142   delete fAttributes;
143   fAttributes=0;
144  }
145 }
146 ///////////////////////////////////////////////////////////////////////////
147 AliCalorimeter::AliCalorimeter(Int_t nrow,Int_t ncol)
148 {
149 // Create a calorimeter module matrix
150  fNrows=nrow;
151  fNcolumns=ncol;
152  fNsignals=0;
153  fModules=0;
154  fNclusters=0;
155  fClusters=0;
156  fAttributes=new TMatrix(nrow,ncol);
157  fGains=new TMatrix(nrow,ncol);
158  fMatrix=new AliCalmodule**[nrow];
159  fPositions=new AliPosition**[nrow];
160  for (Int_t row=0; row<nrow; row++)
161  {
162   fMatrix[row]=new AliCalmodule*[ncol];
163   fPositions[row]=new AliPosition*[ncol];
164   // Initialise the various matrices
165   for (Int_t col=0; col<ncol; col++)
166   {
167    fMatrix[row][col]=0;
168    fPositions[row][col]=0;
169    (*fGains)(row,col)=1;
170    (*fAttributes)(row,col)=0;
171   }
172  }
173
174  // Mark the edge modules
175  for (Int_t j=0; j<ncol; j++)
176  {
177   (*fAttributes)(0,j)=10;
178   (*fAttributes)(nrow-1,j)=10;
179  }
180  for (Int_t i=0; i<nrow; i++)
181  {
182   (*fAttributes)(i,0)=10;
183   (*fAttributes)(i,ncol-1)=10;
184  }
185  
186  fHmodules=0;
187  fHclusters=0;
188
189  fNvetos=0;
190  fVetos=0;
191
192  fName=" ";
193 }
194 ///////////////////////////////////////////////////////////////////////////
195 Int_t AliCalorimeter::GetNrows()
196 {
197 // Provide the number of rows for the calorimeter module matrix
198  return fNrows;
199 }
200 ///////////////////////////////////////////////////////////////////////////
201 Int_t AliCalorimeter::GetNcolumns()
202 {
203 // Provide the number of columns for the calorimeter module matrix
204  return fNcolumns;
205 }
206 ///////////////////////////////////////////////////////////////////////////
207 void AliCalorimeter::SetSignal(Int_t row,Int_t col,Float_t sig)
208 {
209 // Set the signal for a certain calorimeter module
210  
211  if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
212  
213  if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
214  {
215   AliCalmodule* m=fMatrix[row-1][col-1];
216   if (!m) // only count new modules
217   {
218    if (!fModules)
219    {
220     fModules=new TObjArray();  // Default size, expanded automatically
221     fModules->SetOwner();
222    }
223    fNsignals++;
224    m=new AliCalmodule();
225    AliPosition* r=fPositions[row-1][col-1];
226    if (r) m->SetPosition(*r);
227    fModules->Add(m);
228    fMatrix[row-1][col-1]=m;
229   }
230   m->SetSignal(row,col,sig);
231  }
232  else
233  {
234   cout << " *AliCalorimeter::SetSignal* row,col : " << row << "," << col
235        << " out of range." << endl;
236   cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
237  }
238 }
239 ///////////////////////////////////////////////////////////////////////////
240 void AliCalorimeter::AddSignal(Int_t row, Int_t col, Float_t sig)
241 {
242 // Add the signal to a certain calorimeter module
243  
244  if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
245  
246  if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
247  {
248   AliCalmodule* m=fMatrix[row-1][col-1];
249   if (!m) // only count new modules
250   {
251    SetSignal(row,col,sig);
252   }
253   else
254   {
255    m->AddSignal(row,col,sig);
256   }
257  }
258  else
259  {
260   cout << " *AliCalorimeter::AddSignal* row,col : " << row << "," << col
261        << " out of range." << endl;
262   cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
263  }
264 }
265 ///////////////////////////////////////////////////////////////////////////
266 void AliCalorimeter::AddSignal(AliCalmodule* mod)
267 {
268 // Add the signal of module mod to the current calorimeter data.
269 // This enables mixing of calorimeter data of various events.
270  
271  if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
272  
273  Int_t row=mod->GetRow();
274  Int_t col=mod->GetColumn();
275  Float_t sig=mod->GetSignal();
276  AliPosition r=mod->GetPosition();
277
278  if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
279  {
280   AliCalmodule* m=fMatrix[row-1][col-1];
281   if (!m) // No module existed yet at this position
282   {
283    if (!fModules)
284    {
285     fModules=new TObjArray();  // Default size, expanded automatically
286     fModules->SetOwner();
287    }
288    fNsignals++;
289    m=new AliCalmodule;
290    fModules->Add(m);
291    fMatrix[row-1][col-1]=m;
292    m->SetPosition(r);
293   }
294   m->AddSignal(row,col,sig);
295   if (!fPositions[row-1][col-1]) fPositions[row-1][col-1]=new AliPosition(r);
296  }
297  else
298  {
299   cout << " *AliCalorimeter::AddSignal* row,col : " << row << "," << col
300        << " out of range." << endl;
301   cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
302  }
303 }
304 ///////////////////////////////////////////////////////////////////////////
305 void AliCalorimeter::Reset(Int_t row,Int_t col)
306 {
307 // Reset the signal for a certain calorimeter module
308  
309  if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
310  
311  if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
312  {
313   AliCalmodule* m=fMatrix[row-1][col-1];
314   if (m)
315   {
316    fModules->Remove(m);
317    fNsignals--;
318    fModules->Compress();
319    delete m;
320    fMatrix[row-1][col-1]=0;
321   }
322  }
323  else
324  {
325   cout << " *AliCalorimeter::Reset* row,col : " << row << "," << col
326        << " out of range." << endl;
327   cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
328  }
329 }
330 ///////////////////////////////////////////////////////////////////////////
331 void AliCalorimeter::Reset(Int_t mode)
332 {
333 // Reset the signals for the complete calorimeter.
334 // Normally this is done to prepare for the data of the next event.
335 //
336 // mode = 0 : Module positions, gains, edge and dead flags remain unchanged.
337 //        1 : Module positions, gains, edge and dead flags are cleared.
338 //
339 // The default is mode=0.
340 //
341 // Note : In the case of reading AliCalorimeter objects from a data file,
342 //        one has to reset the AliCalorimeter object with mode=1
343 //        (or explicitly delete it) before reading-in the next object
344 //        in order to prevent memory leaks.
345
346  if (mode<0 || mode>1)
347  {
348   cout << " *AliCalorimeter::Reset* Wrong argument. mode = " << mode << endl;
349   return;
350  }
351
352  fNsignals=0;
353  if (fModules)
354  {
355   delete fModules;
356   fModules=0;
357  }
358
359  fNclusters=0;
360  if (fClusters)
361  {
362   delete fClusters;
363   fClusters=0;
364  }
365
366  fNvetos=0;
367  if (fVetos)
368  {
369   delete fVetos;
370   fVetos=0;
371  }
372
373  if (fMatrix || fPositions)
374  {
375   for (Int_t i=0; i<fNrows; i++)
376   {
377    for (Int_t j=0; j<fNcolumns; j++)
378    {
379     if (mode==0)
380     {
381      fMatrix[i][j]=0;
382     }
383     else
384     {
385      if (fPositions[i][j]) delete fPositions[i][j];
386     }
387    }
388    if (mode==1)
389    {
390     if (fPositions[i]) delete [] fPositions[i];
391     if (fMatrix[i]) delete [] fMatrix[i];
392    }
393   }
394  }
395
396  // Free memory allocated for various arrays and matrices.
397  if (mode==1)
398  {
399   if (fMatrix)
400   {
401    delete [] fMatrix;
402    fMatrix=0;
403   }
404   if (fPositions)
405   {
406    delete [] fPositions;
407    fPositions=0;
408   }
409   if (fGains)
410   {
411    delete fGains;
412    fGains=0;
413   }
414   if (fAttributes)
415   {
416    delete fAttributes;
417    fAttributes=0;
418   }
419   if (fHmodules)
420   {
421    delete fHmodules;
422    fHmodules=0;
423   }
424   if (fHclusters)
425   {
426    delete fHclusters;
427    fHclusters=0;
428   }
429  }
430 }
431 ///////////////////////////////////////////////////////////////////////////
432 Float_t AliCalorimeter::GetSignal(Int_t row,Int_t col)
433 {
434 // Provide the signal of a certain calorimeter module.
435 // In case the module was marked dead, 0 is returned.
436  
437  if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
438  
439  if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
440  {
441   Float_t signal=0;
442   AliCalmodule* m=fMatrix[row-1][col-1];
443   if (m)
444   {
445    Int_t dead=m->GetDeadValue();
446    if (!dead) signal=m->GetSignal();
447   }
448   return signal;
449  }
450  else
451  {
452   cout << " *AliCalorimeter::GetSignal* row,col : " << row << "," << col
453        << " out of range." << endl;
454   cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
455   return 0;
456  }
457 }
458 ///////////////////////////////////////////////////////////////////////////
459 void AliCalorimeter::SetEdgeOn(Int_t row,Int_t col)
460 {
461 // Indicate a certain calorimeter module as 'edge module'
462  
463  if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
464  
465  if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
466  {
467   Float_t word=(*fAttributes)(row-1,col-1);
468   Int_t iword=int(word+0.1);
469   Int_t dead=iword%10;
470   Int_t edge=1;
471   (*fAttributes)(row-1,col-1)=float(dead+10*edge); 
472  }
473  else
474  {
475   cout << " *AliCalorimeter::SetEdgeOn* row,col : " << row << "," << col
476        << " out of range." << endl;
477   cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
478  }
479 }
480 ///////////////////////////////////////////////////////////////////////////
481 void AliCalorimeter::SetEdgeOff(Int_t row,Int_t col)
482 {
483 // Indicate a certain calorimeter module as 'non-edge module'
484  
485  if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
486  
487  if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
488  {
489   Float_t word=(*fAttributes)(row-1,col-1);
490   Int_t iword=int(word+0.1);
491   Int_t dead=iword%10;
492   Int_t edge=0;
493   (*fAttributes)(row-1,col-1)=float(dead+10*edge); 
494  }
495  else
496  {
497   cout << " *AliCalorimeter::SetEdgeOff* row,col : " << row << "," << col
498        << " out of range." << endl;
499   cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
500  }
501 }
502 ///////////////////////////////////////////////////////////////////////////
503 void AliCalorimeter::SetDead(Int_t row,Int_t col)
504 {
505 // Indicate a certain calorimeter module as 'dead module'
506  
507  if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
508  
509  if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
510  {
511   Float_t word=(*fAttributes)(row-1,col-1);
512   Int_t iword=int(word+0.1);
513   Int_t edge=iword/10;
514   Int_t dead=1;
515   (*fAttributes)(row-1,col-1)=float(dead+10*edge);
516   if (fMatrix[row-1][col-1]) (fMatrix[row-1][col-1])->SetDead();
517  
518   // Increase the 'edge value' of surrounding modules
519   Int_t rlow=row-1;
520   Int_t rup=row+1;
521   Int_t clow=col-1;
522   Int_t cup=col+1;
523  
524   if (rlow < 1) rlow=row;
525   if (rup > fNrows) rup=fNrows;
526   if (clow < 1) clow=col;
527   if (cup > fNcolumns) cup=fNcolumns;
528  
529   for (Int_t i=rlow; i<=rup; i++)
530   {
531    for (Int_t j=clow; j<=cup; j++)
532    {
533     if (i!=row || j!=col) // No increase of edge value for the dead module itself
534     {
535      word=(*fAttributes)(i-1,j-1);
536      iword=int(word+0.1);
537      edge=iword/10;
538      dead=iword%10;
539      edge++;
540      (*fAttributes)(i-1,j-1)=float(dead+10*edge);
541     }
542    }
543   }
544  }
545  else
546  {
547   cout << " *AliCalorimeter::SetDead* row,col : " << row << "," << col
548        << " out of range." << endl;
549   cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
550  }
551 }
552 ///////////////////////////////////////////////////////////////////////////
553 void AliCalorimeter::SetAlive(Int_t row,Int_t col)
554 {
555 // Indicate a certain calorimeter module as 'active module'
556  
557  if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
558  
559  if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
560  {
561   Float_t word=(*fAttributes)(row-1,col-1);
562   Int_t iword=int(word+0.1);
563   Int_t edge=iword/10;
564   Int_t dead=0;
565   (*fAttributes)(row-1,col-1)=float(dead+10*edge);
566   if (fMatrix[row-1][col-1]) (fMatrix[row-1][col-1])->SetAlive();
567  
568   // Decrease the 'edge value' of surrounding modules
569   Int_t rlow=row-1;
570   Int_t rup=row+1;
571   Int_t clow=col-1;
572   Int_t cup=col+1;
573  
574   if (rlow < 1) rlow=row;
575   if (rup > fNrows) rup=fNrows;
576   if (clow < 1) clow=col;
577   if (cup > fNcolumns) cup=fNcolumns;
578  
579   for (Int_t i=rlow; i<=rup; i++)
580   {
581    for (Int_t j=clow; j<=cup; j++)
582    {
583     if (i!=row || j!=col) // No decrease of edge value for the dead module itself
584     {
585      word=(*fAttributes)(i-1,j-1);
586      iword=int(word+0.1);
587      edge=iword/10;
588      dead=iword%10;
589      if (edge>0) edge--;
590      (*fAttributes)(i-1,j-1)=float(dead+10*edge);
591     }
592    }
593   }
594  }
595  else
596  {
597   cout << " *AliCalorimeter::SetAlive* row,col : " << row << "," << col
598        << " out of range." << endl;
599   cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
600  }
601 }
602 ///////////////////////////////////////////////////////////////////////////
603 void AliCalorimeter::SetGain(Int_t row,Int_t col,Float_t gain)
604 {
605 // Set the gain value for a certain calorimeter module
606  
607  if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
608  
609  if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
610  {
611   (*fGains)(row-1,col-1)=gain;
612   if (fMatrix[row-1][col-1]) (fMatrix[row-1][col-1])->SetGain(gain);
613  }
614  else
615  {
616   cout << " *AliCalorimeter::SetGain* row,col : " << row << "," << col
617        << " out of range." << endl;
618   cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
619  }
620 }
621 ///////////////////////////////////////////////////////////////////////////
622 void AliCalorimeter::SetPosition(Int_t row,Int_t col,Float_t* vec,TString f)
623 {
624 // Set the position in user coordinates for a certain calorimeter module
625  
626  if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
627  
628  if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
629  {
630   if (!fPositions[row-1][col-1]) fPositions[row-1][col-1]=new AliPosition;
631   (fPositions[row-1][col-1])->SetVector(vec,f);
632   if (fMatrix[row-1][col-1]) (fMatrix[row-1][col-1])->SetPosition(vec,f);
633  }
634  else
635  {
636   cout << " *AliCalorimeter::SetPosition* row,col : " << row << "," << col
637        << " out of range." << endl;
638   cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
639  }
640 }
641 ///////////////////////////////////////////////////////////////////////////
642 void AliCalorimeter::SetPosition(Int_t row,Int_t col,Ali3Vector& r)
643 {
644 // Set the position for a certain calorimeter module
645  
646  if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
647  
648  if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
649  {
650   if (!fPositions[row-1][col-1]) fPositions[row-1][col-1]=new AliPosition;
651   (fPositions[row-1][col-1])->SetPosition(r);
652   if (fMatrix[row-1][col-1]) (fMatrix[row-1][col-1])->SetPosition(r);
653  }
654  else
655  {
656   cout << " *AliCalorimeter::SetPosition* row,col : " << row << "," << col
657        << " out of range." << endl;
658   cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
659  }
660 }
661 ///////////////////////////////////////////////////////////////////////////
662 Int_t AliCalorimeter::GetEdgeValue(Int_t row,Int_t col)
663 {
664 // Provide the value of the edge flag of a certain module
665  
666  if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
667  {
668   Float_t word=(*fAttributes)(row-1,col-1);
669   Int_t iword=int(word+0.1);
670   Int_t edge=iword/10;
671   return edge;
672  }
673  else
674  {
675   cout << " *AliCalorimeter::GetEdgeValue* row,col : " << row << "," << col
676        << " out of range." << endl;
677   cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
678   return 0;
679  }
680 }
681 ///////////////////////////////////////////////////////////////////////////
682 Int_t AliCalorimeter::GetDeadValue(Int_t row,Int_t col)
683 {
684 // Provide the value of the dead flag of a certain module
685  
686  if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
687  {
688   Float_t word=(*fAttributes)(row-1,col-1);
689   Int_t iword=int(word+0.1);
690   Int_t dead=iword%10;
691   return dead;
692  }
693  else
694  {
695   cout << " *AliCalorimeter::GetDeadValue* row,col : " << row << "," << col
696        << " out of range." << endl;
697   cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
698   return 0;
699  }
700 }
701 ///////////////////////////////////////////////////////////////////////////
702 Float_t AliCalorimeter::GetGain(Int_t row,Int_t col)
703 {
704 // Provide the gain value of a certain module
705  
706  if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
707  {
708    return (*fGains)(row-1,col-1);
709  }
710  else
711  {
712   cout << " *AliCalorimeter::GetGain* row,col : " << row << "," << col
713        << " out of range." << endl;
714   cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
715   return 0;
716  }
717 }
718 ///////////////////////////////////////////////////////////////////////////
719 void AliCalorimeter::GetPosition(Int_t row,Int_t col,Float_t* vec,TString f)
720 {
721 // Return the position in user coordinates for a certain calorimeter module
722  
723  if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
724  
725  if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
726  {
727 //  if (fMatrix[row-1][col-1]) (fMatrix[row-1][col-1])->GetPosition(vec,f);
728   if (fPositions[row-1][col-1]) (fPositions[row-1][col-1])->GetVector(vec,f);
729  }
730  else
731  {
732   cout << " *AliCalorimeter::GetPosition* row,col : " << row << "," << col
733        << " out of range." << endl;
734   cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
735  }
736 }
737 ///////////////////////////////////////////////////////////////////////////
738 AliPosition* AliCalorimeter::GetPosition(Int_t row,Int_t col)
739 {
740 // Access to the position of a certain calorimeter module
741  
742  if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
743  
744  if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
745  {
746   return fPositions[row-1][col-1];
747  }
748  else
749  {
750   cout << " *AliCalorimeter::GetPosition* row,col : " << row << "," << col
751        << " out of range." << endl;
752   cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
753   return 0;
754  }
755 }
756 ///////////////////////////////////////////////////////////////////////////
757 Float_t AliCalorimeter::GetClusteredSignal(Int_t row,Int_t col)
758 {
759 // Provide the module signal after clustering
760  
761  if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
762  
763  if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
764  {
765   if (fMatrix[row-1][col-1])
766   {
767    return (fMatrix[row-1][col-1])->GetClusteredSignal();
768   }
769   else
770   {
771    return 0;
772   }
773  }
774  else
775  {
776   cout << " *AliCalorimeter::GetClusteredSignal* row,col : " << row << "," << col
777        << " out of range." << endl;
778   cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
779   return 0;
780  }
781 }
782 ///////////////////////////////////////////////////////////////////////////
783 Int_t AliCalorimeter::GetNsignals()
784 {
785 // Provide the number of modules that contain a signal
786 // Note : The number of modules marked 'dead' but which had a signal
787 //        are included.
788  return fNsignals;
789 }
790 ///////////////////////////////////////////////////////////////////////////
791 void AliCalorimeter::Group(Int_t n)
792 {
793 // Group the individual modules into clusters
794 // Module signals of n rings around the central module will be grouped
795  
796  if (fNsignals > 0) // Directly return if no modules fired
797  {
798   if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
799  
800   if (fNclusters > 0) Ungroup(); // Restore unclustered situation if needed
801  
802   // Order the modules with decreasing signal
803   AliCalmodule** ordered=new AliCalmodule*[fNsignals]; // temp. array for ordered modules
804   Int_t nord=0;
805   Sortm(ordered,nord);
806  
807   // Clustering of modules. Start with the highest signal.
808   if (fClusters)
809   {
810    delete fClusters;
811    fClusters=0;
812   }
813   fClusters=new TObjArray();
814   fClusters->SetOwner();
815   fNclusters=0;
816   Int_t row=0;
817   Int_t col=0;
818   AliCalcluster* c=0;
819   AliCalmodule* m=0;
820   for (Int_t i=0; i<nord; i++)
821   {
822    row=ordered[i]->GetRow();    // row number of cluster center
823    col=ordered[i]->GetColumn(); // column number of cluster center
824    if (row>0 && row<=fNrows && col>0 && col<=fNcolumns)
825    {
826     m=fMatrix[row-1][col-1];
827     if (!m) continue;
828
829     // only use modules not yet used in a cluster
830     if (m->GetClusteredSignal() > 0.)
831     {
832      Int_t edge=GetEdgeValue(row,col);
833      c=new AliCalcluster;
834      if (!edge) c->Start(*m);   // module to start the cluster if not on edge
835      if (c->GetNmodules() > 0)  // cluster started successfully (no edge)
836      {
837       fClusters->Add(c);
838       fNclusters++;       // update cluster counter
839       AddRing(row,col,n); // add signals of n rings around the center
840      }
841      else
842      {
843       if (c) delete c;
844       c=0;
845      }
846     }
847    }
848   }
849  
850   // Delete the temp. array
851   if (ordered)
852   { 
853    for (Int_t j=0; j<nord; j++)
854    {
855     ordered[j]=0;
856    }
857    delete [] ordered;
858   }
859  }
860 }
861 ///////////////////////////////////////////////////////////////////////////
862 void AliCalorimeter::Sortm(AliCalmodule** ordered,Int_t& nord)
863 {
864 // Order the modules with decreasing signal
865  
866  nord=0;
867  for (Int_t i=0; i<fNrows; i++) // loop over all modules of the matrix
868  {
869   for (Int_t ii=0; ii<fNcolumns; ii++)
870   {
871    if (GetSignal(i+1,ii+1) <= 0.) continue; // only take alive modules with a signal
872  
873    if (nord == 0) // store the first module with a signal at the first ordered position
874    {
875     nord++;
876     ordered[nord-1]=fMatrix[i][ii];
877     continue;
878    }
879  
880    for (Int_t j=0; j<=nord; j++) // put module in the right ordered position
881    {
882     if (j == nord) // module has smallest signal seen so far
883     {
884      nord++;
885      ordered[j]=fMatrix[i][ii]; // add module at the end
886      break; // go for next matrix module
887     }
888  
889     if (GetSignal(i+1,ii+1) < ordered[j]->GetSignal()) continue;
890  
891     nord++;
892     for (Int_t k=nord-1; k>j; k--) // create empty position
893     {
894      ordered[k]=ordered[k-1];
895     }
896     ordered[j]=fMatrix[i][ii]; // put module at empty position
897     break; // go for next matrix module
898    }
899   }
900  }
901 }
902 ///////////////////////////////////////////////////////////////////////////
903 void AliCalorimeter::AddRing(Int_t row, Int_t col, Int_t n)
904 {
905 // Add module signals of 1 ring around (row,col) to current cluster
906 // n denotes the maximum number of rings around cluster center
907 // Note : This function is used recursively
908  
909  if (n >= 1) // Check if any rings left for recursive calls
910  {
911   Float_t signal=GetSignal(row,col); // signal of (row,col) module
912  
913   Int_t lrow=row-1; if (lrow < 1) lrow=1;                 // row lowerbound for ring
914   Int_t urow=row+1; if (urow > fNrows) urow=fNrows;       // row upperbound for ring
915   Int_t lcol=col-1; if (lcol < 1) lcol=1;                 // col lowerbound for ring
916   Int_t ucol=col+1; if (ucol > fNcolumns) ucol=fNcolumns; // row upperbound for ring
917  
918   for (Int_t i=lrow; i<=urow; i++)
919   {
920    for (Int_t j=lcol; j<=ucol; j++)
921    {
922     // add module(i,j) to cluster if the signal <= signal(row,col)
923     if (GetSignal(i,j) <= signal)
924     {
925      AliCalmodule* m=fMatrix[i-1][j-1];
926      if (m) ((AliCalcluster*)fClusters->At(fNclusters-1))->Add(*m);
927     }
928     AddRing(i,j,n-1); // Go for ring of modules around this (i,j) one
929    }
930   }
931  }
932 }
933 ///////////////////////////////////////////////////////////////////////////
934 Int_t AliCalorimeter::GetNclusters()
935 {
936 // Provide the number of clusters
937  return fNclusters;
938 }
939 ///////////////////////////////////////////////////////////////////////////
940 AliCalcluster* AliCalorimeter::GetCluster(Int_t j)
941 {
942 // Provide cluster number j
943 // Note : j=1 denotes the first cluster
944  if ((j >= 1) && (j <= fNclusters))
945  {
946   return (AliCalcluster*)fClusters->At(j-1);
947  }
948  else
949  {
950   cout << " *AliCalorimeter::GetCluster* cluster number : " << j
951        << " out of range." << endl;
952   cout << " -- Cluster number 1 (if any) returned " << endl;
953   return (AliCalcluster*)fClusters->At(0);
954  }
955 }
956 ///////////////////////////////////////////////////////////////////////////
957 AliCalmodule* AliCalorimeter::GetModule(Int_t j)
958 {
959 // Provide 'fired' module number j
960 // Note : j=1 denotes the first 'fired' module
961  if ((j >= 1) && (j <= fNsignals))
962  {
963   return (AliCalmodule*)fModules->At(j-1);
964  }
965  else
966  {
967   cout << " *AliCalorimeter::GetModule* module number : " << j
968        << " out of range." << endl;
969   cout << " -- Fired module number 1 (if any) returned " << endl;
970   return (AliCalmodule*)fModules->At(0);
971  }
972 }
973 ///////////////////////////////////////////////////////////////////////////
974 AliCalmodule* AliCalorimeter::GetModule(Int_t row,Int_t col)
975 {
976 // Provide access to module (row,col).
977 // Note : first module is at (1,1).
978
979  if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
980
981  if (row>=1 && row<=fNrows && col>=1 && col<=fNcolumns)
982  {
983   return fMatrix[row-1][col-1];
984  }
985  else
986  {
987   cout << " *AliCalorimeter::GetModule* row,col : " << row << ", " << col
988        << " out of range." << endl;
989   return 0;
990  }
991 }
992 ///////////////////////////////////////////////////////////////////////////
993 TH2F* AliCalorimeter::DrawModules()
994 {
995 // Provide a lego plot of the module signals
996  
997  if (fHmodules)
998  {
999   fHmodules->Reset();
1000  }
1001  else
1002  {
1003   fHmodules=new TH2F("fHmodules","Module signals",
1004             fNcolumns,0.5,float(fNcolumns)+0.5,fNrows,0.5,float(fNrows)+0.5);
1005  
1006   fHmodules->SetDirectory(0); // Suppress global character of histo pointer
1007  }
1008  
1009  AliCalmodule* m;
1010  Float_t row,col,signal;
1011  Int_t dead;
1012  for (Int_t i=0; i<fNsignals; i++)
1013  {
1014   m=(AliCalmodule*)fModules->At(i);
1015   if (m)
1016   {
1017    row=float(m->GetRow());
1018    col=float(m->GetColumn());
1019    dead=m->GetDeadValue();
1020    signal=0;
1021    if (!dead) signal=m->GetSignal();
1022    if (signal>0.) fHmodules->Fill(col,row,signal);
1023   }
1024  }
1025  
1026  fHmodules->Draw("lego");
1027  return fHmodules;
1028 }
1029 ///////////////////////////////////////////////////////////////////////////
1030 TH2F* AliCalorimeter::DrawClusters()
1031 {
1032 // Provide a lego plot of the cluster signals
1033  
1034  if (fHclusters)
1035  {
1036   fHclusters->Reset();
1037  }
1038  else
1039  {
1040   fHclusters=new TH2F("fHclusters","Cluster signals",
1041             fNcolumns,0.5,float(fNcolumns)+0.5,fNrows,0.5,float(fNrows)+0.5);
1042  
1043   fHclusters->SetDirectory(0); // Suppress global character of histo pointer
1044  }
1045  
1046  AliCalcluster* c;
1047  Float_t row,col,signal;
1048  for (Int_t i=0; i<fNclusters; i++)
1049  {
1050   c=(AliCalcluster*)fClusters->At(i);
1051   if (c)
1052   {
1053    row=float(c->GetRow());
1054    col=float(c->GetColumn());
1055    signal=c->GetSignal();
1056    if (signal>0.) fHclusters->Fill(col,row,signal);
1057   }
1058  }
1059  
1060  fHclusters->Draw("lego");
1061  return fHclusters;
1062 }
1063 ///////////////////////////////////////////////////////////////////////////
1064 void AliCalorimeter::LoadMatrix()
1065 {
1066 // Load the Calorimeter module matrix data back from the TObjArray
1067
1068  // Initialise the module matrix
1069  if (!fMatrix)
1070  {
1071   fMatrix=new AliCalmodule**[fNrows];
1072   for (Int_t i=0; i<fNrows; i++)
1073   {
1074    fMatrix[i]=new AliCalmodule*[fNcolumns];
1075   }
1076  }
1077
1078  // Initialise the position matrix
1079  if (!fPositions)
1080  {
1081   fPositions=new AliPosition**[fNrows];
1082   for (Int_t j=0; j<fNrows; j++)
1083   {
1084    fPositions[j]=new AliPosition*[fNcolumns];
1085   }
1086  }
1087
1088  for (Int_t jrow=0; jrow<fNrows; jrow++)
1089  {
1090   for (Int_t jcol=0; jcol<fNcolumns; jcol++)
1091   {
1092    fMatrix[jrow][jcol]=0;
1093    fPositions[jrow][jcol]=0;
1094   }
1095  }
1096  
1097  // Copy the module pointers back into the matrix
1098  AliCalmodule* m=0;
1099  Int_t row=0;
1100  Int_t col=0;
1101  Int_t nsig=0;
1102  if (fModules) nsig=fModules->GetEntries();
1103  for (Int_t j=0; j<nsig; j++)
1104  {
1105   m=(AliCalmodule*)fModules->At(j);
1106   if (m)
1107   {
1108    row=m->GetRow();
1109    col=m->GetColumn();
1110    AliPosition r=m->GetPosition();
1111    fMatrix[row-1][col-1]=m;
1112    fPositions[row-1][col-1]=new AliPosition(r);
1113   }
1114  }
1115 }
1116 ///////////////////////////////////////////////////////////////////////////
1117 void AliCalorimeter::Ungroup()
1118 {
1119 // Set the module signals back to the non-clustered situation
1120  
1121  if (!fMatrix) LoadMatrix(); // Restore matrix data in case of reading input
1122  
1123  Int_t nsig=0;
1124  if (fModules) nsig=fModules->GetEntries();
1125
1126  Float_t signal=0;
1127  AliCalmodule* m=0;
1128  for (Int_t j=0; j<nsig; j++)
1129  {
1130   m=(AliCalmodule*)fModules->At(j);
1131   if (m)
1132   {
1133    signal=m->GetSignal();
1134    m->SetClusteredSignal(signal);
1135   }
1136  }
1137 }
1138 ///////////////////////////////////////////////////////////////////////////
1139 void AliCalorimeter::AddVetoSignal(AliSignal& s)
1140 {
1141 // Associate an (extrapolated) AliSignal as veto to the calorimeter.
1142  if (!fVetos)
1143  {
1144   fNvetos=0;
1145   fVetos=new TObjArray();
1146   fVetos->SetOwner();
1147  } 
1148
1149  Int_t nvalues=s.GetNvalues();
1150  AliSignal* sx=new AliSignal(nvalues);
1151  sx->SetName(s.GetName());
1152  
1153  sx->SetPosition((Ali3Vector&)s);
1154
1155  Double_t sig,err;
1156  for (Int_t i=1; i<=nvalues; i++)
1157  {
1158   sig=s.GetSignal(i);
1159   err=s.GetSignalError(i);
1160   sx->SetSignal(sig,i);
1161   sx->SetSignalError(err,i);
1162  } 
1163
1164  fVetos->Add(sx);
1165  fNvetos++;
1166 }
1167 ///////////////////////////////////////////////////////////////////////////
1168 Int_t AliCalorimeter::GetNvetos()
1169 {
1170 // Provide the number of veto signals associated to the calorimeter
1171  return fNvetos;
1172 }
1173 ///////////////////////////////////////////////////////////////////////////
1174 AliSignal* AliCalorimeter::GetVetoSignal(Int_t i)
1175 {
1176 // Provide access to the i-th veto signal of this calorimeter
1177 // Note : The first hit corresponds to i=1
1178
1179  if (i>0 && i<=fNvetos)
1180  {
1181   return (AliSignal*)fVetos->At(i-1);
1182  }
1183  else
1184  {
1185   cout << " *AliCalorimeter::GetVetoSignal* Signal number " << i
1186        << " out of range." << endl;
1187   cout << " --- First signal (if any) returned." << endl;
1188   return (AliSignal*)fVetos->At(0);
1189  }
1190 }
1191 ///////////////////////////////////////////////////////////////////////////
1192 void AliCalorimeter::SetName(TString name)
1193 {
1194 // Set the name of the calorimeter system.
1195  fName=name;
1196 }
1197 ///////////////////////////////////////////////////////////////////////////
1198 TString AliCalorimeter::GetName()
1199 {
1200 // Provide the name of the calorimeter system.
1201  return fName;
1202 }
1203 ///////////////////////////////////////////////////////////////////////////