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