]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RALICE/AliCalorimeter.cxx
03-oct-2003 NvE Typos fixed in AliCalorimeter.cxx.
[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 generic 2D geometry is used in which a module is identified by (row,col).
22 // Obviously this geometry can be a matrix, but also any other regular
23 // structure is supported, provided the user has adopted a proper convention
24 // to uniquely address a module via the (row,col) indices.  
25 // Note : First module is identified as (1,1).
26 //
27 // This is the way to define and enter signals into a calorimeter :
28 //
29 //   AliCalorimeter cal;
30 //
31 //   cal.AddSignal(5,7,85.4);
32 //   cal.AddSignal(5,7,25.9);
33 //   cal.AddSignal(3,5,1000);
34 //   cal.SetSignal(5,7,10.3);
35 //   cal.Reset(3,5);             // Reset module (3,5) as being 'not fired'
36 //                               // All module data are re-initialised.
37 //   cal.SetEdgeOn(1,1);         // Declare module (1,1) as an 'edge module'
38 //   cal.SetDead(8,3);
39 //   cal.SetGain(2,8,3.2);
40 //
41 //   Float_t vec[3]={6,1,20};
42 //   cal.SetPosition(2,8,vec,"car");
43 //
44 //   AliSignal s;
45 //   Float_t loc[3]={-1,12,3};
46 //   s.SetPosition(loc,"car");
47 //   s.SetSignal(328);
48 //   cal.AddVetoSignal(s); // Associate (extrapolated) signal as a veto
49 //
50 //   cal.Group(2);      // Group 'fired' modules into clusters
51 //                      // Perform grouping over 2 rings around the center
52 //   cal.Reset();       // Reset the complete calorimeter
53 //                      // Normally to prepare for the next event data
54 //                      // Note : Module gain, offset, edge and dead flags remain
55 //
56 //--- Author: Nick van Eijndhoven 13-jun-1997 UU-SAP Utrecht
57 //- Modified: NvE $Date$ UU-SAP Utrecht
58 ///////////////////////////////////////////////////////////////////////////
59
60 #include "AliCalorimeter.h"
61 #include "Riostream.h"
62
63 ClassImp(AliCalorimeter) // Class implementation to enable ROOT I/O
64  
65 AliCalorimeter::AliCalorimeter() : TObject()
66 {
67 // Default constructor, all parameters set to 0.
68 // Create a calorimeter module matrix with fixed row and column size.
69 // Note : Due to the dynamic size extension when signals are set,
70 //        the  "edge modules" can NOT be marked automatically.
71 //        This has to be done manually by the user via the SetEdgeOn() 
72 //        memberfunction.
73  fNrows=0;
74  fNcolumns=0;
75  fSwap=0;
76  fMatrix=0;
77  fClusters=0;
78  fHmodules=0;
79  fHclusters=0;
80  fVetos=0;
81  fAttributes=0;
82  fPositions=0;
83  fName="Unspecified AliCalorimeter";
84 }
85 ///////////////////////////////////////////////////////////////////////////
86 AliCalorimeter::~AliCalorimeter()
87 {
88 // Destructor to delete memory allocated to the various arrays and matrices
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  if (fMatrix)
110  {
111   delete fMatrix;
112   fMatrix=0;
113  }
114  if (fPositions)
115  {
116   delete fPositions;
117   fPositions=0;
118  }
119  if (fAttributes)
120  {
121   delete fAttributes;
122   fAttributes=0;
123  }
124 }
125 ///////////////////////////////////////////////////////////////////////////
126 AliCalorimeter::AliCalorimeter(Int_t nrow,Int_t ncol) : TObject()
127 {
128 // Create a calorimeter module matrix with fixed row and column size.
129 // The modules at the edges are automatically marked as "edge modules".
130  fNrows=nrow;
131  fNcolumns=ncol;
132  fClusters=0;
133
134  fSwap=0;
135  fMatrix=0;
136  fPositions=0;
137
138  fAttributes=new TObjArray(nrow);
139  fAttributes->SetOwner();
140  
141  // Mark the edge modules
142  for (Int_t row=1; row<=nrow; row++)
143  {
144   AliAttribObj* a=new AliAttribObj();
145   if (row==1 || row==nrow)
146   {
147    for (Int_t col=1; col<=ncol; col++)
148    {
149     a->SetEdgeOn(col);
150    }
151   }
152   else
153   {
154    a->SetEdgeOn(1);
155    a->SetEdgeOn(ncol);
156   }
157   fAttributes->Add(a);
158  }
159  
160  fHmodules=0;
161  fHclusters=0;
162
163  fVetos=0;
164
165  fName="Unspecified AliCalorimeter";
166 }
167 ///////////////////////////////////////////////////////////////////////////
168 AliCalorimeter::AliCalorimeter(AliCalorimeter& c) : TObject(c)
169 {
170 // Copy constructor
171  fClusters=0;
172  fVetos=0;
173
174  fAttributes=0;
175
176  fHmodules=0;
177  fHclusters=0;
178
179  fMatrix=0;
180  fPositions=0;
181
182  fNrows=c.fNrows;
183  fNcolumns=c.fNcolumns;
184  fName=c.fName;
185
186  fSwap=c.fSwap;
187
188  if (c.fPositions)
189  {
190   Int_t nrows=(c.fPositions)->GetMaxRow();
191   Int_t ncols=(c.fPositions)->GetMaxColumn();
192   for (Int_t irow=1; irow<=nrows; irow++)
193   {
194    for (Int_t icol=1; icol<=ncols; icol++)
195    {
196     AliPosition* p=c.GetPosition(irow,icol);
197     if (p) SetPosition(irow,icol,*p);
198    } 
199   }  
200  }
201
202  Int_t size=0;
203  Int_t n=0;
204  if (c.fAttributes)
205  {
206   size=c.fAttributes->GetSize();
207   n=c.fAttributes->GetEntries();
208  }
209  if (size)
210  {
211   fAttributes=new TObjArray(size);
212   fAttributes->SetOwner();
213   for (Int_t ia=0; ia<n; ia++)
214   {
215    AliAttribObj* a=(AliAttribObj*)(c.fAttributes->At(ia));
216    if (a) fAttributes->Add(new AliAttribObj(*a));
217   }
218  }
219
220  n=c.GetNclusters();
221  if (n)
222  {
223   fClusters=new TObjArray();
224   fClusters->SetOwner();
225   for (Int_t icl=1; icl<=n; icl++)
226   {
227    AliCalcluster* cl=c.GetCluster(icl);
228    if (cl) fClusters->Add(new AliCalcluster(*cl));
229   }
230  }
231
232  n=c.GetNsignals();
233  if (n)
234  {
235   fMatrix=new AliObjMatrix();
236   fMatrix->SetOwner();
237   fMatrix->SetSwapMode(fSwap);
238   for (Int_t im=1; im<=n; im++)
239   {
240    AliCalmodule* m=c.GetModule(im);
241    if (m) fMatrix->EnterObject(m->GetRow(),m->GetColumn(),new AliCalmodule(*m));
242   }
243  }
244
245  n=c.GetNvetos();
246  for (Int_t iv=1; iv<=n; iv++)
247  {
248   AliSignal* s=c.GetVetoSignal(iv);
249   if (s) AddVetoSignal(s);
250  }
251 }
252 ///////////////////////////////////////////////////////////////////////////
253 Int_t AliCalorimeter::GetNrows()
254 {
255 // Provide the number of rows for the calorimeter module matrix
256  Int_t nrows=fNrows;
257  if (fMatrix && !nrows) nrows=fMatrix->GetMaxRow();
258  return nrows;
259 }
260 ///////////////////////////////////////////////////////////////////////////
261 Int_t AliCalorimeter::GetNcolumns()
262 {
263 // Provide the number of columns for the calorimeter module matrix
264  Int_t ncols=fNcolumns;
265  if (fMatrix && !ncols) ncols=fMatrix->GetMaxColumn();
266  return ncols;
267 }
268 ///////////////////////////////////////////////////////////////////////////
269 void AliCalorimeter::SetSignal(Int_t row,Int_t col,Float_t sig)
270 {
271 // Set the signal for a certain calorimeter module.
272
273  // Check for (row,col) boundaries.
274  if (row<1 || col<1 || (fNrows && fNcolumns && (row>fNrows || col>fNcolumns)))
275  {
276   cout << " *AliCalorimeter::SetSignal* row,col : " << row << "," << col
277        << " out of range." << endl;
278   if (fNrows && fNcolumns) cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
279   return;
280  }
281  
282  if (!fMatrix)
283  {
284   fMatrix=new AliObjMatrix();
285   fMatrix->SetOwner();
286   fMatrix->SetSwapMode(fSwap);
287  }
288
289  AliCalmodule* m=GetModule(row,col);
290  if (!m) // initialise for a new module
291  {
292   m=new AliCalmodule();
293   AliPosition* r=0;
294   if (fPositions) r=(AliPositionObj*)fPositions->GetObject(row,col);
295   if (r) m->SetPosition(*r);
296   if (fAttributes)
297   {
298    AliAttribObj* a=0;
299    if (row <= fAttributes->GetSize()) a=(AliAttribObj*)fAttributes->At(row-1);
300    if (a)
301    {
302     if (a->GetGainFlag(col)) m->SetGain(a->GetGain(col));
303     if (a->GetOffsetFlag(col)) m->SetOffset(a->GetOffset(col));
304     if (a->GetDeadValue(col)) m->SetDead();
305     if (a->GetEdgeValue(col)) m->SetEdgeValue(a->GetEdgeValue(col));
306    }
307   }
308   fMatrix->EnterObject(row,col,m);
309  }
310
311  m->SetSignal(row,col,sig);
312 }
313 ///////////////////////////////////////////////////////////////////////////
314 void AliCalorimeter::AddSignal(Int_t row, Int_t col, Float_t sig)
315 {
316 // Add the signal to a certain calorimeter module.
317
318  // Check for (row,col) boundaries.
319  if (row<1 || col<1 || (fNrows && fNcolumns && (row>fNrows || col>fNcolumns)))
320  {
321   cout << " *AliCalorimeter::AddSignal* row,col : " << row << "," << col
322        << " out of range." << endl;
323   if (fNrows && fNcolumns) cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
324   return;
325  }
326  
327   AliCalmodule* m=GetModule(row,col);
328   if (!m) // initialise for new modules
329   {
330    SetSignal(row,col,sig);
331   }
332   else
333   {
334    m->AddSignal(row,col,sig);
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 // Note : The position and attributes according to the user provided data
344 //        for the corresponding (row,col) location will be used.
345 //        In case there is no user provided data present, the position and
346 //        attributes of the first module added to the corresponding (row,col)
347 //        location will be taken, except for the "edge" and "dead" indicators.
348 //        The latter will then both be set to 0.
349
350  if (!mod) return;
351  
352  Int_t row=mod->GetRow();
353  Int_t col=mod->GetColumn();
354  Float_t sig=mod->GetSignal();
355
356  // Check for (row,col) boundaries.
357  if (row<1 || col<1 || (fNrows && fNcolumns && (row>fNrows || col>fNcolumns)))
358  {
359   cout << " *AliCalorimeter::AddSignal* row,col : " << row << "," << col
360        << " out of range." << endl;
361   if (fNrows && fNcolumns) cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
362   return;
363  }
364  
365  if (!fMatrix)
366  {
367   fMatrix=new AliObjMatrix();
368   fMatrix->SetOwner();
369   fMatrix->SetSwapMode(fSwap);
370  }
371
372  AliCalmodule* m=GetModule(row,col);
373  if (!m) // No module existed yet at this position
374  {
375   m=new AliCalmodule(*mod);
376   AliPosition* r=0;
377   if (fPositions) r=(AliPositionObj*)fPositions->GetObject(row,col);
378   if (r) m->SetPosition(*r);
379   // Don't take the dead and edge attributes from this module,
380   // but from the calorimeter dbase, if present.
381   m->SetEdgeOff();
382   m->SetAlive();
383   if (fAttributes)
384   {
385    AliAttribObj* a=0;
386    if (row <= fAttributes->GetSize()) a=(AliAttribObj*)fAttributes->At(row-1);
387    if (a)
388    {
389     if (a->GetGainFlag(col)) m->SetGain(a->GetGain(col));
390     if (a->GetOffsetFlag(col)) m->SetOffset(a->GetOffset(col));
391     if (a->GetDeadValue(col)) m->SetDead();
392     if (a->GetEdgeValue(col)) m->SetEdgeValue(a->GetEdgeValue(col));
393    }
394   }
395   fMatrix->EnterObject(row,col,m);
396  }
397  else
398  {
399   m->AddSignal(row,col,sig);
400  }
401 }
402 ///////////////////////////////////////////////////////////////////////////
403 void AliCalorimeter::Reset(Int_t row,Int_t col)
404 {
405 // Reset the signal for a certain calorimeter module.
406 // Note : Module position and attributes remain unchanged.
407
408  // Check for (row,col) boundaries.
409  if (row<1 || col<1 || (fNrows && fNcolumns && (row>fNrows || col>fNcolumns)))
410  {
411   cout << " *AliCalorimeter::Reset* row,col : " << row << "," << col
412        << " out of range." << endl;
413   if (fNrows && fNcolumns) cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
414   return;
415  }
416  
417  AliCalmodule* m=GetModule(row,col);
418  if (m) fMatrix->RemoveObject(row,col);
419 }
420 ///////////////////////////////////////////////////////////////////////////
421 void AliCalorimeter::Reset(Int_t mode)
422 {
423 // Reset the signals for the complete calorimeter.
424 // Normally this is done to prepare for the data of the next event.
425 //
426 // mode = 0 : Swap mode, module positions and attributes remain unchanged.
427 //        1 : Swap mode, module positions and attributes are cleared.
428 //
429 // The default is mode=0.
430 //
431 // Note : In the case of reading AliCalorimeter objects from a data file,
432 //        one has to reset the AliCalorimeter object with mode=1
433 //        (or explicitly delete it) before reading-in the next object
434 //        in order to prevent memory leaks.
435
436  if (mode<0 || mode>1)
437  {
438   cout << " *AliCalorimeter::Reset* Wrong argument. mode = " << mode << endl;
439   return;
440  }
441
442  if (fClusters)
443  {
444   delete fClusters;
445   fClusters=0;
446  }
447
448  if (fVetos)
449  {
450   delete fVetos;
451   fVetos=0;
452  }
453
454  if (mode==1)
455  {
456   if (fMatrix)
457   {
458    delete fMatrix;
459    fMatrix=0;
460   }
461   if (fPositions)
462   {
463    delete fPositions;
464    fPositions=0;
465   }
466  }
467  else
468  {
469   if (fMatrix) fMatrix->Reset();
470  }
471
472  // Free memory allocated for the various arrays.
473  if (mode==1)
474  {
475   if (fAttributes)
476   {
477    delete fAttributes;
478    fAttributes=0;
479   }
480   if (fHmodules)
481   {
482    delete fHmodules;
483    fHmodules=0;
484   }
485   if (fHclusters)
486   {
487    delete fHclusters;
488    fHclusters=0;
489   }
490  }
491 }
492 ///////////////////////////////////////////////////////////////////////////
493 Float_t AliCalorimeter::GetSignal(Int_t row,Int_t col,Int_t mode)
494 {
495 // Provide the signal of a certain calorimeter module.
496 // In case the module was marked dead, 0 is returned.
497 //
498 // mode = 0 : Just the module signal is returned
499 //        1 : The module signal is corrected for the gain and offset.
500 //            In case the gain value was not set, gain=1 will be assumed.
501 //            In case the gain value was 0, a signal value of 0 is returned.
502 //            In case the offset value was not set, offset=0 will be assumed.
503 //
504 // The corrected signal (sigc) is determined as follows :
505 //
506 //              sigc=(signal/gain)-offset 
507 //
508 // The default is mode=0.
509
510  // Check for (row,col) boundaries.
511  if (row<1 || col<1 || (fNrows && fNcolumns && (row>fNrows || col>fNcolumns)))
512  {
513   cout << " *AliCalorimeter::GetSignal* row,col : " << row << "," << col
514        << " out of range." << endl;
515   if (fNrows && fNcolumns) cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
516   return 0;
517  }
518  
519  Float_t signal=0;
520  Float_t gain=1;
521  Float_t offset=0;
522  AliCalmodule* m=GetModule(row,col);
523  if (m)
524  {
525   Int_t dead=m->GetDeadValue();
526   if (!dead) signal=m->GetSignal();
527
528   if (mode==0 || dead) return signal;
529
530   // Correct the signal for the gain and offset
531   if (GetGainFlag(row,col))
532   {
533    gain=GetGain(row,col);
534   }
535   else
536   {
537    if (m->GetGainFlag()) gain=m->GetGain();
538   }
539
540   if (GetOffsetFlag(row,col))
541   {
542    offset=GetOffset(row,col);
543   }
544   else
545   {
546    if (m->GetOffsetFlag()) offset=m->GetOffset();
547   }
548
549   if (fabs(gain)>0.)
550   {
551    signal=(signal/gain)-offset;
552   }
553   else
554   {
555    signal=0;
556   }
557  }
558  return signal;
559 }
560 ///////////////////////////////////////////////////////////////////////////
561 void AliCalorimeter::SetEdgeOn(Int_t row,Int_t col)
562 {
563 // Indicate a certain calorimeter module as 'edge module'.
564
565  // Check for (row,col) boundaries.
566  if (row<1 || col<1 || (fNrows && fNcolumns && (row>fNrows || col>fNcolumns)))
567  {
568   cout << " *AliCalorimeter::SetEdgeOn* row,col : " << row << "," << col
569        << " out of range." << endl;
570   if (fNrows && fNcolumns) cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
571   return;
572  }
573
574  if (!fAttributes)
575  {
576   fAttributes=new TObjArray(row);
577   fAttributes->SetOwner();
578  }
579  else
580  {
581   if (row > fAttributes->GetSize()) fAttributes->Expand(row);
582  }
583
584  AliAttribObj* a=(AliAttribObj*)fAttributes->At(row-1);
585  if (a)
586  {
587   a->SetEdgeOn(col);
588  }
589  else
590  {
591   a=new AliAttribObj();
592   a->SetEdgeOn(col);
593   fAttributes->AddAt(a,row-1);
594  } 
595
596  AliCalmodule* m=GetModule(row,col);
597  if (m) m->SetEdgeOn();
598 }
599 ///////////////////////////////////////////////////////////////////////////
600 void AliCalorimeter::SetEdgeOff(Int_t row,Int_t col)
601 {
602 // Indicate a certain calorimeter module as 'non-edge module'.
603
604  // Check for (row,col) boundaries.
605  if (row<1 || col<1 || (fNrows && fNcolumns && (row>fNrows || col>fNcolumns)))
606  {
607   cout << " *AliCalorimeter::SetEdgeOff* row,col : " << row << "," << col
608        << " out of range." << endl;
609   if (fNrows && fNcolumns) cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
610   return;
611  }
612
613  // Only action on fAttributes in case an attribute is present at (row,col),
614  // since by default a module has edge=0 unless explicitly set otherwise.
615  if (fAttributes)
616  {
617   if (row <= fAttributes->GetSize())
618   {
619    AliAttribObj* a=(AliAttribObj*)fAttributes->At(row-1);
620    if (a) a->SetEdgeOff(col);
621   }
622  } 
623
624  AliCalmodule* m=GetModule(row,col);
625  if (m) m->SetEdgeOff();
626 }
627 ///////////////////////////////////////////////////////////////////////////
628 void AliCalorimeter::SetDead(Int_t row,Int_t col)
629 {
630 // Indicate a certain calorimeter module as 'dead module'
631
632  // Check for (row,col) boundaries in case of a fixed size calorimeter
633  if (row<1 || col<1 || (fNrows && fNcolumns && (row>fNrows || col>fNcolumns)))
634  {
635   cout << " *AliCalorimeter::SetDead* row,col : " << row << "," << col
636        << " out of range." << endl;
637   if (fNrows && fNcolumns) cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
638   return;
639  }
640
641  // Make Attributes storage 1 row (and also 1 column) larger than needed
642  // because the 'edge value' of the (future) surrounding modules has
643  // to be updated as well.
644  if (!fAttributes)
645  {
646   fAttributes=new TObjArray(row+1);
647   fAttributes->SetOwner();
648  }
649  else
650  {
651   if (row >= fAttributes->GetSize()) fAttributes->Expand(row+1);
652  }
653
654  AliAttribObj* a=(AliAttribObj*)fAttributes->At(row-1);
655  if (a)
656  {
657   a->SetDead(col);
658  }
659  else
660  {
661   a=new AliAttribObj();
662   a->SetDead(col);
663   fAttributes->AddAt(a,row-1);
664  } 
665
666  AliCalmodule* m=GetModule(row,col);
667  if (m) m->SetDead();
668  
669  // Increase the 'edge value' of surrounding modules
670  Int_t rlow=row-1;
671  Int_t rup=row+1;
672  Int_t clow=col-1;
673  Int_t cup=col+1;
674  
675  if (rlow < 1) rlow=row;
676  if (clow < 1) clow=col;
677  
678  for (Int_t i=rlow; i<=rup; i++)
679  {
680   for (Int_t j=clow; j<=cup; j++)
681   {
682    if (i!=row || j!=col) // No increase of edge value for the 'dead' module itself
683    {
684     AliAttribObj* a=(AliAttribObj*)fAttributes->At(i-1);
685     if (a)
686     {
687      a->IncreaseEdgeValue(j);
688     }
689     else
690     {
691     a=new AliAttribObj();
692     a->SetEdgeOn(j);
693     fAttributes->AddAt(a,i-1);
694     } 
695
696     AliCalmodule* m=GetModule(i,j);
697     if (m) m->IncreaseEdgeValue();
698    }
699   }
700  }
701 }
702 ///////////////////////////////////////////////////////////////////////////
703 void AliCalorimeter::SetAlive(Int_t row,Int_t col)
704 {
705 // Indicate a certain calorimeter module as 'active module'.
706
707  // Check for (row,col) boundaries.
708  if (row<1 || col<1 || (fNrows && fNcolumns && (row>fNrows || col>fNcolumns)))
709  {
710   cout << " *AliCalorimeter::SetAlive* row,col : " << row << "," << col
711        << " out of range." << endl;
712   if (fNrows && fNcolumns) cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
713   return;
714  }
715
716  // Only action on fAttributes in case an attribute is present at (row,col),
717  // since by default a module has dead=0 unless explicitly set otherwise.
718  if (fAttributes)
719  {
720   if (row <= fAttributes->GetSize())
721   {
722    AliAttribObj* a=(AliAttribObj*)fAttributes->At(row-1);
723    if (a) a->SetAlive(col);
724   }
725  } 
726
727  AliCalmodule* m=GetModule(row,col);
728  if (m) m->SetAlive();
729  
730  // Decrease the 'edge value' of surrounding modules
731  Int_t rlow=row-1;
732  Int_t rup=row+1;
733  Int_t clow=col-1;
734  Int_t cup=col+1;
735  
736  if (rlow < 1) rlow=row;
737  if (clow < 1) clow=col;
738
739  for (Int_t i=rlow; i<=rup; i++)
740  {
741   for (Int_t j=clow; j<=cup; j++)
742   {
743    if (i!=row || j!=col) // No decrease of edge value for the 'alive' module itself
744    {
745     if (i <= fAttributes->GetSize())
746     {
747      AliAttribObj* a=(AliAttribObj*)fAttributes->At(i-1);
748      if (a) a->DecreaseEdgeValue(j);
749     }
750     AliCalmodule* m=GetModule(i,j);
751     if (m) m->DecreaseEdgeValue();
752    }
753   }
754  }
755 }
756 ///////////////////////////////////////////////////////////////////////////
757 void AliCalorimeter::SetGain(Int_t row,Int_t col,Float_t gain)
758 {
759 // Set the gain value for a certain calorimeter module.
760 // See the memberfunction GetSignal() for a definition of the gain value.
761
762  // Check for (row,col) boundaries.
763  if (row<1 || col<1 || (fNrows && fNcolumns && (row>fNrows || col>fNcolumns)))
764  {
765   cout << " *AliCalorimeter::SetGain* row,col : " << row << "," << col
766        << " out of range." << endl;
767   if (fNrows && fNcolumns) cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
768   return;
769  }
770
771  if (!fAttributes)
772  {
773   fAttributes=new TObjArray(row);
774   fAttributes->SetOwner();
775  }
776  else
777  {
778   if (row > fAttributes->GetSize()) fAttributes->Expand(row);
779  }
780
781  AliAttribObj* a=(AliAttribObj*)fAttributes->At(row-1);
782  if (a)
783  {
784   a->SetGain(gain,col);
785  }
786  else
787  {
788   a=new AliAttribObj();
789   a->SetGain(gain,col);
790   fAttributes->AddAt(a,row-1);
791  } 
792
793  AliCalmodule* m=GetModule(row,col);
794  if (m) m->SetGain(gain);
795 }
796 ///////////////////////////////////////////////////////////////////////////
797 void AliCalorimeter::SetOffset(Int_t row,Int_t col,Float_t offset)
798 {
799 // Set the offset value for a certain calorimeter module.
800 // See the memberfunction GetSignal() for a definition of the offset value.
801
802  // Check for (row,col) boundaries.
803  if (row<1 || col<1 || (fNrows && fNcolumns && (row>fNrows || col>fNcolumns)))
804  {
805   cout << " *AliCalorimeter::SetOffset* row,col : " << row << "," << col
806        << " out of range." << endl;
807   if (fNrows && fNcolumns) cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
808   return;
809  }
810
811  if (!fAttributes)
812  {
813   fAttributes=new TObjArray(row);
814   fAttributes->SetOwner();
815  }
816  else
817  {
818   if (row > fAttributes->GetSize()) fAttributes->Expand(row);
819  }
820
821  AliAttribObj* a=(AliAttribObj*)fAttributes->At(row-1);
822  if (a)
823  {
824   a->SetOffset(offset,col);
825  }
826  else
827  {
828   a=new AliAttribObj();
829   a->SetOffset(offset,col);
830   fAttributes->AddAt(a,row-1);
831  } 
832
833  AliCalmodule* m=GetModule(row,col);
834  if (m) m->SetOffset(offset);
835 }
836 ///////////////////////////////////////////////////////////////////////////
837 void AliCalorimeter::SetPosition(Int_t row,Int_t col,Float_t* vec,TString f)
838 {
839 // Set the position in user coordinates for a certain calorimeter module
840  Ali3Vector r;
841  r.SetVector(vec,f);
842  SetPosition(row,col,r);
843
844 ///////////////////////////////////////////////////////////////////////////
845 void AliCalorimeter::SetPosition(Int_t row,Int_t col,Ali3Vector& r)
846 {
847 // Set the position for a certain calorimeter module
848
849  // Check for (row,col) boundaries.
850  if (row<1 || col<1 || (fNrows && fNcolumns && (row>fNrows || col>fNcolumns)))
851  {
852   cout << " *AliCalorimeter::SetPosition* row,col : " << row << "," << col
853        << " out of range." << endl;
854   if (fNrows && fNcolumns) cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
855   return;
856  }
857
858  if (!fPositions)
859  {
860   fPositions=new AliObjMatrix();
861   fPositions->SetOwner();
862   fPositions->SetSwapMode(fSwap);
863  }
864
865  AliPositionObj* p=(AliPositionObj*)fPositions->GetObject(row,col);
866
867  if (p)
868  {
869   p->Load(r);
870  }
871  else
872  {
873   p=new AliPositionObj();
874   p->Load(r);
875   fPositions->EnterObject(row,col,p);
876  }
877
878  // Update the position of the calorimeter module itself as well if it exists
879  AliCalmodule* m=GetModule(row,col);
880  if (m) m->SetPosition(r);
881 }
882 ///////////////////////////////////////////////////////////////////////////
883 Int_t AliCalorimeter::GetEdgeValue(Int_t row,Int_t col)
884 {
885 // Provide the value of the edge flag of a certain module.
886
887  // Check for (row,col) boundaries.
888  if (row<1 || col<1 || (fNrows && fNcolumns && (row>fNrows || col>fNcolumns)))
889  {
890   cout << " *AliCalorimeter::GetEdgeValue* row,col : " << row << "," << col
891        << " out of range." << endl;
892   if (fNrows && fNcolumns) cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
893   return 0;
894  }
895
896  Int_t edge=0;
897
898  if (fAttributes)
899  {
900   if (row <= fAttributes->GetSize())
901   {
902    AliAttribObj* a=(AliAttribObj*)fAttributes->At(row-1);
903    if (a)
904    {
905     if (col <= a->GetNcalflags())
906     {
907      edge=a->GetEdgeValue(col);
908      return edge;
909     }
910    }
911   }
912  }
913
914  AliCalmodule* m=GetModule(row,col);
915  if (m) edge=m->GetEdgeValue();
916  return edge;
917 }
918 ///////////////////////////////////////////////////////////////////////////
919 Int_t AliCalorimeter::GetDeadValue(Int_t row,Int_t col)
920 {
921 // Provide the value of the dead flag of a certain module
922
923  // Check for (row,col) boundaries.
924  if (row<1 || col<1 || (fNrows && fNcolumns && (row>fNrows || col>fNcolumns)))
925  {
926   cout << " *AliCalorimeter::GetDeadValue* row,col : " << row << "," << col
927        << " out of range." << endl;
928   if (fNrows && fNcolumns) cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
929   return 0;
930  }
931
932  Int_t dead=0;
933
934  if (fAttributes)
935  {
936   if (row <= fAttributes->GetSize())
937   {
938    AliAttribObj* a=(AliAttribObj*)fAttributes->At(row-1);
939    if (a)
940    {
941     if (col <= a->GetNcalflags())
942     {
943      dead=a->GetDeadValue(col);
944      return dead;
945     }
946    }
947   }
948  }
949
950  AliCalmodule* m=GetModule(row,col);
951  if (m) dead=m->GetDeadValue();
952  return dead;
953 }
954 ///////////////////////////////////////////////////////////////////////////
955 Int_t AliCalorimeter::GetGainFlag(Int_t row,Int_t col)
956 {
957 // Provide the value of the gain flag of a certain module.
958
959  // Check for (row,col) boundaries.
960  if (row<1 || col<1 || (fNrows && fNcolumns && (row>fNrows || col>fNcolumns)))
961  {
962   cout << " *AliCalorimeter::GetGainFlag* row,col : " << row << "," << col
963        << " out of range." << endl;
964   if (fNrows && fNcolumns) cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
965   return 0;
966  }
967
968  Int_t gf=0;
969
970  if (fAttributes)
971  {
972   if (row <= fAttributes->GetSize())
973   {
974    AliAttribObj* a=(AliAttribObj*)fAttributes->At(row-1);
975    if (a)
976    {
977     if (col <= a->GetNcalflags())
978     {
979      gf=a->GetGainFlag(col);
980      return gf;
981     }
982    }
983   }
984  }
985
986  AliCalmodule* m=GetModule(row,col);
987  if (m) gf=m->GetGainFlag();
988  return gf;
989 }
990 ///////////////////////////////////////////////////////////////////////////
991 Int_t AliCalorimeter::GetOffsetFlag(Int_t row,Int_t col)
992 {
993 // Provide the value of the offset flag of a certain module.
994
995  // Check for (row,col) boundaries.
996  if (row<1 || col<1 || (fNrows && fNcolumns && (row>fNrows || col>fNcolumns)))
997  {
998   cout << " *AliCalorimeter::GetOffsetFlag* row,col : " << row << "," << col
999        << " out of range." << endl;
1000   if (fNrows && fNcolumns) cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
1001   return 0;
1002  }
1003
1004  Int_t of=0;
1005
1006  if (fAttributes)
1007  {
1008   if (row <= fAttributes->GetSize())
1009   {
1010    AliAttribObj* a=(AliAttribObj*)fAttributes->At(row-1);
1011    if (a)
1012    {
1013     if (col <= a->GetNcalflags())
1014     {
1015      of=a->GetOffsetFlag(col);
1016      return of;
1017     }
1018    }
1019   }
1020  }
1021
1022  AliCalmodule* m=GetModule(row,col);
1023  if (m) of=m->GetOffsetFlag();
1024  return of;
1025 }
1026 ///////////////////////////////////////////////////////////////////////////
1027 Float_t AliCalorimeter::GetGain(Int_t row,Int_t col)
1028 {
1029 // Provide the gain value of a certain module.
1030 // See the memberfunction GetSignal() for a definition of the gain value.
1031 //
1032 // In case the gain value is unknown, the value 0 will be returned. 
1033
1034  // Check for (row,col) boundaries.
1035  if (row<1 || col<1 || (fNrows && fNcolumns && (row>fNrows || col>fNcolumns)))
1036  {
1037   cout << " *AliCalorimeter::GetGain* row,col : " << row << "," << col
1038        << " out of range." << endl;
1039   if (fNrows && fNcolumns) cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
1040   return 0;
1041  }
1042
1043  Float_t gain=0;
1044
1045  if (fAttributes)
1046  {
1047   if (row <= fAttributes->GetSize())
1048   {
1049    AliAttribObj* a=(AliAttribObj*)fAttributes->At(row-1);
1050    if (a)
1051    {
1052     if (col <= a->GetNcalflags())
1053     {
1054      if (a->GetGainFlag(col))
1055      {
1056       gain=a->GetGain(col);
1057       return gain;
1058      }
1059     }
1060    }
1061   }
1062  }
1063
1064  AliCalmodule* m=GetModule(row,col);
1065  if (m)
1066  {
1067   if (m->GetGainFlag())
1068   {
1069    gain=m->GetGain();
1070   }
1071  }
1072  return gain;
1073 }
1074 ///////////////////////////////////////////////////////////////////////////
1075 Float_t AliCalorimeter::GetOffset(Int_t row,Int_t col)
1076 {
1077 // Provide the offset value of a certain module.
1078 // See the memberfunction GetSignal() for a definition of the offset value.
1079 //
1080 // In case the offset value is unknown, the value 0 will be returned. 
1081
1082  // Check for (row,col) boundaries.
1083  if (row<1 || col<1 || (fNrows && fNcolumns && (row>fNrows || col>fNcolumns)))
1084  {
1085   cout << " *AliCalorimeter::GetOffset* row,col : " << row << "," << col
1086        << " out of range." << endl;
1087   if (fNrows && fNcolumns) cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
1088   return 0;
1089  }
1090
1091  Float_t offset=0;
1092
1093  if (fAttributes)
1094  {
1095   if (row <= fAttributes->GetSize())
1096   {
1097    AliAttribObj* a=(AliAttribObj*)fAttributes->At(row-1);
1098    if (a)
1099    {
1100     if (col <= a->GetNcalflags())
1101     {
1102      if (a->GetOffsetFlag(col))
1103      {
1104       offset=a->GetOffset(col);
1105       return offset;
1106      }
1107     }
1108    }
1109   }
1110  }
1111
1112  AliCalmodule* m=GetModule(row,col);
1113  if (m)
1114  {
1115   if (m->GetOffsetFlag())
1116   {
1117    offset=m->GetOffset();
1118   }
1119  }
1120  return offset;
1121 }
1122 ///////////////////////////////////////////////////////////////////////////
1123 void AliCalorimeter::GetPosition(Int_t row,Int_t col,Float_t* vec,TString f)
1124 {
1125 // Return the position in user coordinates for a certain calorimeter module
1126  vec[0]=0;
1127  vec[1]=0;
1128  vec[2]=0;
1129
1130  AliPosition* p=GetPosition(row,col);
1131  if (p) p->GetVector(vec,f);
1132 }
1133 ///////////////////////////////////////////////////////////////////////////
1134 AliPosition* AliCalorimeter::GetPosition(Int_t row,Int_t col)
1135 {
1136 // Access to the position of a certain calorimeter module.
1137
1138  // Check for (row,col) boundaries.
1139  if (row<1 || col<1 || (fNrows && fNcolumns && (row>fNrows || col>fNcolumns)))
1140  {
1141   cout << " *AliCalorimeter::GetPosition* row,col : " << row << "," << col
1142        << " out of range." << endl;
1143   if (fNrows && fNcolumns) cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
1144   return 0;
1145  }
1146
1147  if (!fPositions && !fMatrix) return 0;
1148
1149  AliPositionObj* po=0;
1150  if (fPositions) po=(AliPositionObj*)fPositions->GetObject(row,col);
1151  if (po) return po;
1152
1153  AliCalmodule* m=GetModule(row,col);
1154  return m;
1155 }
1156 ///////////////////////////////////////////////////////////////////////////
1157 Float_t AliCalorimeter::GetClusteredSignal(Int_t row,Int_t col)
1158 {
1159 // Provide the module signal after clustering.
1160
1161  // Check for (row,col) boundaries.
1162  if (row<1 || col<1 || (fNrows && fNcolumns && (row>fNrows || col>fNcolumns)))
1163  {
1164   cout << " *AliCalorimeter::GetClusteredSignal* row,col : " << row << "," << col
1165        << " out of range." << endl;
1166   if (fNrows && fNcolumns) cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
1167   return 0;
1168  }
1169  
1170  Float_t sig=0;
1171
1172  AliCalmodule* m=GetModule(row,col);
1173  if (m) sig=m->GetClusteredSignal();
1174
1175  return sig;
1176 }
1177 ///////////////////////////////////////////////////////////////////////////
1178 Int_t AliCalorimeter::GetNsignals()
1179 {
1180 // Provide the number of modules that contain a signal
1181 // Note : The number of modules marked 'dead' but which had a signal
1182 //        are included.
1183  Int_t nsig=0;
1184  if (fMatrix) nsig=fMatrix->GetNobjects();
1185  return nsig;
1186 }
1187 ///////////////////////////////////////////////////////////////////////////
1188 void AliCalorimeter::Group(Int_t n,Int_t mode)
1189 {
1190 // Group the individual modules into clusters.
1191 // Module signals of n rings around the central module will be grouped.
1192 // The grouping process will start with the module containing the highest signal
1193 // in an iterative way.
1194 // For this all fired modules are ordered w.r.t. decreasing signal.
1195 // The search mode for the module signal hierarchy can be specified by the user.
1196 //
1197 // mode = 1 : Search performed via the (row,col) structure of the matrix (SortM)
1198 //        2 : Search performed via the linear array of fired modules (SortA)
1199 //
1200 // See the docs of the memberfunctions SortM and SortA for additional details.
1201 //
1202 // Default values : n=1 mode=1.
1203
1204  if (mode<1 || mode>2)
1205  {
1206   cout << " *AliCalorimeter::Group* Invalid mode : " << mode << endl;
1207   cout << " Default value mode=1 will be used." << endl;
1208   mode=1;
1209  }
1210
1211  if (!fMatrix) return;
1212  
1213  Int_t nsignals=GetNsignals();
1214  if (nsignals > 0) // Only continue if there are fired modules
1215  {
1216   if (GetNclusters() > 0) Ungroup(); // Restore unclustered situation if needed
1217  
1218   // Order the modules with decreasing signal
1219   AliCalmodule** ordered=new AliCalmodule*[nsignals]; // temp. array for ordered modules
1220   Int_t nord=0;
1221   if (mode==1) SortM(ordered,nord);
1222   if (mode==2) SortA(ordered,nord);
1223  
1224   // Clustering of modules. Start with the highest signal.
1225   if (fClusters)
1226   {
1227    delete fClusters;
1228    fClusters=0;
1229   }
1230   fClusters=new TObjArray();
1231   fClusters->SetOwner();
1232   Int_t row=0;
1233   Int_t col=0;
1234   AliCalcluster* c=0;
1235   AliCalmodule* m=0;
1236   for (Int_t i=0; i<nord; i++)
1237   {
1238    row=ordered[i]->GetRow();    // row number of cluster center
1239    col=ordered[i]->GetColumn(); // column number of cluster center
1240
1241    m=(AliCalmodule*)fMatrix->GetObject(row,col);
1242    if (!m) continue;
1243
1244    // only use modules not yet used in a cluster
1245    if (m->GetClusteredSignal() > 0.)
1246    {
1247     Int_t edge=GetEdgeValue(row,col);
1248     c=new AliCalcluster();
1249     if (!edge) c->Start(*m);   // module to start the cluster if not on edge
1250     if (c->GetNmodules() > 0)  // cluster started successfully (no edge)
1251     {
1252      fClusters->Add(c);
1253      AddRing(row,col,n); // add signals of n rings around the center
1254     }
1255     else
1256     {
1257      if (c) delete c;
1258      c=0;
1259     }
1260    }
1261   }
1262  
1263   // Delete the temp. array
1264   if (ordered)
1265   { 
1266    for (Int_t j=0; j<nord; j++)
1267    {
1268     ordered[j]=0;
1269    }
1270    delete [] ordered;
1271   }
1272  }
1273 }
1274 ///////////////////////////////////////////////////////////////////////////
1275 void AliCalorimeter::SortM(AliCalmodule** ordered,Int_t& nord)
1276 {
1277 // Order the modules with decreasing signal by looping over the (row,col) grid
1278 // of the matrix.
1279 // Modules which were declared as "Dead" will be rejected.
1280 // The gain etc... corrected module signals will be used in the ordering process.
1281 //
1282 // Note : This method may become slow for large, very finely granulated calorimeters.
1283 //
1284 // Very specific case :
1285 // ====================
1286 // In case of various overlapping showers of which the central modules have
1287 // EXACTLY the same signal this ordering procedure may have the following
1288 // advantages and disadvantages.
1289 //
1290 // Advantages :
1291 // ------------
1292 // * In case of multi-overlapping showers, the central shower will NOT
1293 //   be "eaten-up" from both sides, resulting in a slightly more accurate
1294 //   cluster signal.
1295 // * This method produces re-producable results, irrespective of the filling
1296 //   order of the matrix modules.
1297 //
1298 // Disadvantages :
1299 // ---------------
1300 // * In case of a very high occupancy, there might be a slight effect on the
1301 //   cluster signals depending on the geometrical location in the detector matrix.
1302  
1303  Int_t nrows=fMatrix->GetMaxRow();
1304  Int_t ncols=fMatrix->GetMaxColumn();
1305
1306  Float_t signal=0.;
1307  nord=0;
1308  for (Int_t irow=1; irow<=nrows; irow++) // loop over all modules of the matrix
1309  {
1310   for (Int_t icol=1; icol<=ncols; icol++)
1311   {
1312    signal=GetSignal(irow,icol,1); // get the gain etc... corrected signal
1313    if (signal <= 0.) continue;    // only take alive modules with a signal
1314  
1315    if (nord == 0) // store the first module with a signal at the first ordered position
1316    {
1317     nord++;
1318     ordered[nord-1]=(AliCalmodule*)fMatrix->GetObject(irow,icol);
1319     continue;
1320    }
1321  
1322    for (Int_t j=0; j<=nord; j++) // put module in the right ordered position
1323    {
1324     if (j == nord) // module has smallest signal seen so far
1325     {
1326      nord++;
1327      ordered[j]=(AliCalmodule*)fMatrix->GetObject(irow,icol); // add module at the end
1328      break; // go for next matrix module
1329     }
1330  
1331     if (signal < ordered[j]->GetSignal(1,1)) continue;
1332  
1333     nord++;
1334     for (Int_t k=nord-1; k>j; k--) // create empty position
1335     {
1336      ordered[k]=ordered[k-1];
1337     }
1338     ordered[j]=(AliCalmodule*)fMatrix->GetObject(irow,icol); // put module at empty position
1339     break; // go for next matrix module
1340    }
1341   }
1342  }
1343 }
1344 ///////////////////////////////////////////////////////////////////////////
1345 void AliCalorimeter::SortA(AliCalmodule** ordered,Int_t& nord)
1346 {
1347 // Order the modules with decreasing signal by looping over the linear array
1348 // of fired modules.
1349 // Modules which were declared as "Dead" will be rejected.
1350 // The gain etc... corrected module signals will be used in the ordering process.
1351 //
1352 // Note : This method is rather fast even for large, very finely granulated calorimeters.
1353 //
1354 // Very specific case :
1355 // ====================
1356 // In case of various overlapping showers of which the central modules have
1357 // EXACTLY the same signal this ordering procedure may have the following
1358 // advantages and disadvantages.
1359 //
1360 // Advantages :
1361 // ------------
1362 // * Even in case of a very high occupancy, the resulting cluster signals
1363 //   will in general NOT depend on the geometrical location in the detector matrix.
1364 //
1365 // Disadvantages :
1366 // ---------------
1367 // * In case of multi-overlapping showers, the central shower might be
1368 //   "eaten-up" from both sides, resulting in a slightly too low value
1369 //   of the resulting cluster signal.
1370 // * This method might produce results depending on the filling
1371 //   order of the matrix modules.
1372  
1373  Int_t nmods=0;
1374  if (fMatrix) nmods=fMatrix->GetNobjects();
1375
1376  nord=0;
1377  for (Int_t i=1; i<=nmods; i++) // loop over all fired modules of the matrix
1378  {
1379   AliCalmodule* m=(AliCalmodule*)fMatrix->GetObject(i);
1380
1381   if (!m) continue;
1382   if (m->GetDeadValue()) continue; // only take alive modules with a signal
1383  
1384   if (nord == 0) // store the first module with a signal at the first ordered position
1385   {
1386    nord++;
1387    ordered[nord-1]=m;
1388    continue;
1389   }
1390  
1391   for (Int_t j=0; j<=nord; j++) // put module in the right ordered position
1392   {
1393    if (j == nord) // module has smallest signal seen so far
1394    {
1395     nord++;
1396     ordered[j]=m; // add module at the end
1397     break; // go for next matrix module
1398    }
1399  
1400    if (m->GetSignal(1,1) < ordered[j]->GetSignal(1,1)) continue;
1401  
1402    nord++;
1403    for (Int_t k=nord-1; k>j; k--) // create empty position
1404    {
1405     ordered[k]=ordered[k-1];
1406    }
1407    ordered[j]=m; // put module at empty position
1408    break; // go for next matrix module
1409   }
1410  }
1411 }
1412 ///////////////////////////////////////////////////////////////////////////
1413 void AliCalorimeter::AddRing(Int_t row, Int_t col, Int_t n)
1414 {
1415 // Add module signals of 1 ring around (row,col) to current cluster.
1416 // The gain etc... corrected module signals will be used in this process.
1417 // The parameter n denotes the maximum number of rings around cluster center.
1418 // Note : This function is used recursively.
1419
1420  if (!fMatrix) return;
1421
1422  Int_t nrows=fMatrix->GetMaxRow();
1423  Int_t ncols=fMatrix->GetMaxColumn();
1424  
1425  if (n >= 1) // Check if any rings left for recursive calls
1426  {
1427   Float_t signal=GetSignal(row,col,1); // Gain etc... corrected signal of (row,col) module
1428  
1429   Int_t lrow=row-1; if (lrow < 1) lrow=1;         // row lowerbound for ring
1430   Int_t urow=row+1; if (urow > nrows) urow=nrows; // row upperbound for ring
1431   Int_t lcol=col-1; if (lcol < 1) lcol=1;         // col lowerbound for ring
1432   Int_t ucol=col+1; if (ucol > ncols) ucol=ncols; // row upperbound for ring
1433  
1434   for (Int_t i=lrow; i<=urow; i++)
1435   {
1436    for (Int_t j=lcol; j<=ucol; j++)
1437    {
1438     // add module(i,j) to cluster if the signal <= signal(row,col)
1439     if (GetSignal(i,j,1) <= signal)
1440     {
1441      AliCalmodule* m=(AliCalmodule*)fMatrix->GetObject(i,j);
1442      if (m) ((AliCalcluster*)fClusters->At(GetNclusters()-1))->Add(*m);
1443     }
1444     AddRing(i,j,n-1); // Go for ring of modules around this (i,j) one
1445    }
1446   }
1447  }
1448 }
1449 ///////////////////////////////////////////////////////////////////////////
1450 Int_t AliCalorimeter::GetNclusters()
1451 {
1452 // Provide the number of clusters
1453  Int_t nclu=0;
1454  if (fClusters) nclu=fClusters->GetEntries();
1455  return nclu;
1456 }
1457 ///////////////////////////////////////////////////////////////////////////
1458 AliCalcluster* AliCalorimeter::GetCluster(Int_t j)
1459 {
1460 // Provide cluster number j
1461 // Note : j=1 denotes the first cluster
1462
1463  if (!fClusters) return 0;
1464
1465  if ((j >= 1) && (j <= GetNclusters()))
1466  {
1467   return (AliCalcluster*)fClusters->At(j-1);
1468  }
1469  else
1470  {
1471   cout << " *AliCalorimeter::GetCluster* cluster number : " << j
1472        << " out of range ==> 0 returned." << endl;
1473   return 0;
1474  }
1475 }
1476 ///////////////////////////////////////////////////////////////////////////
1477 AliCalmodule* AliCalorimeter::GetModule(Int_t j)
1478 {
1479 // Provide 'fired' module number j
1480 // Note : j=1 denotes the first 'fired' module
1481
1482  if (!fMatrix) return 0;
1483
1484  if ((j >= 1) && (j <= GetNsignals()))
1485  {
1486   return (AliCalmodule*)fMatrix->GetObject(j);
1487  }
1488  else
1489  {
1490   cout << " *AliCalorimeter::GetModule* module number : " << j
1491        << " out of range ==> 0 returned." << endl;
1492   return 0;
1493  }
1494 }
1495 ///////////////////////////////////////////////////////////////////////////
1496 AliCalmodule* AliCalorimeter::GetModule(Int_t row,Int_t col)
1497 {
1498 // Provide access to module (row,col).
1499 // Note : first module is at (1,1).
1500
1501  AliCalmodule* m=0;
1502  if (fMatrix) m=(AliCalmodule*)fMatrix->GetObject(row,col);
1503  return m;
1504 }
1505 ///////////////////////////////////////////////////////////////////////////
1506 TH2F* AliCalorimeter::DrawModules(Float_t thresh,Int_t mode)
1507 {
1508 // Provide a lego plot of the module signals.
1509 // The input parameter mode (default mode=0) has the same meaning as
1510 // specified in the memberfunction GetSignal(row,col,mode).
1511 // Only modules with a (corrected) signal value above the threshold
1512 // (default thresh=0) will be displayed.
1513
1514  Int_t nrows=fNrows;
1515  Int_t ncols=fNcolumns;
1516
1517  if (fMatrix && !nrows && !ncols)
1518  {
1519   nrows=fMatrix->GetMaxRow();
1520   ncols=fMatrix->GetMaxColumn();
1521  }
1522  
1523  if (fHmodules)
1524  {
1525   fHmodules->Reset();
1526  }
1527  else
1528  {
1529   fHmodules=new TH2F("fHmodules","Module signals",
1530             ncols,0.5,float(ncols)+0.5,nrows,0.5,float(nrows)+0.5);
1531  
1532   fHmodules->SetDirectory(0); // Suppress global character of histo pointer
1533  }
1534  
1535  Int_t nmods=GetNsignals();
1536
1537  Float_t row,col,signal;
1538  Int_t dead;
1539  for (Int_t i=1; i<=nmods; i++)
1540  {
1541   AliCalmodule* m=(AliCalmodule*)fMatrix->GetObject(i);
1542   if (m)
1543   {
1544    row=float(m->GetRow());
1545    col=float(m->GetColumn());
1546    dead=m->GetDeadValue();
1547    signal=0;
1548    if (!dead) signal=GetSignal(row,col,mode);
1549    if (signal>thresh) fHmodules->Fill(col,row,signal);
1550   }
1551  }
1552  
1553  fHmodules->Draw("lego");
1554  return fHmodules;
1555 }
1556 ///////////////////////////////////////////////////////////////////////////
1557 TH2F* AliCalorimeter::DrawClusters(Float_t thresh)
1558 {
1559 // Provide a lego plot of the cluster signals.
1560 // Only clusters with a signal value above the threshold (default thresh=0)
1561 // will be displayed.
1562
1563  Int_t nrows=fNrows;
1564  Int_t ncols=fNcolumns;
1565
1566  if (fMatrix && !nrows && !ncols)
1567  {
1568   nrows=fMatrix->GetMaxRow();
1569   ncols=fMatrix->GetMaxColumn();
1570  }
1571  
1572  if (fHclusters)
1573  {
1574   fHclusters->Reset();
1575  }
1576  else
1577  {
1578   fHclusters=new TH2F("fHclusters","Cluster signals",
1579             ncols,0.5,float(ncols)+0.5,nrows,0.5,float(nrows)+0.5);
1580  
1581   fHclusters->SetDirectory(0); // Suppress global character of histo pointer
1582  }
1583  
1584  AliCalcluster* c;
1585  Float_t row,col,signal;
1586  for (Int_t i=0; i<GetNclusters(); i++)
1587  {
1588   c=(AliCalcluster*)fClusters->At(i);
1589   if (c)
1590   {
1591    row=float(c->GetRow());
1592    col=float(c->GetColumn());
1593    signal=c->GetSignal();
1594    if (signal>thresh) fHclusters->Fill(col,row,signal);
1595   }
1596  }
1597  
1598  fHclusters->Draw("lego");
1599  return fHclusters;
1600 }
1601 ///////////////////////////////////////////////////////////////////////////
1602 void AliCalorimeter::Ungroup()
1603 {
1604 // Set the module signals back to the non-clustered situation
1605
1606  if (!fMatrix) return;
1607  
1608  Int_t nsig=GetNsignals();
1609
1610  Float_t signal=0;
1611  for (Int_t j=1; j<=nsig; j++)
1612  {
1613   AliCalmodule* m=(AliCalmodule*)fMatrix->GetObject(j);
1614   if (m)
1615   {
1616    signal=m->GetSignal();
1617    m->SetClusteredSignal(signal);
1618   }
1619  }
1620 }
1621 ///////////////////////////////////////////////////////////////////////////
1622 void AliCalorimeter::AddVetoSignal(AliSignal& s)
1623 {
1624 // Associate an (extrapolated) AliSignal as veto to the calorimeter.
1625  if (!fVetos)
1626  {
1627   fVetos=new TObjArray();
1628   fVetos->SetOwner();
1629  } 
1630
1631  AliSignal* sx=new AliSignal(s);
1632
1633  fVetos->Add(sx);
1634 }
1635 ///////////////////////////////////////////////////////////////////////////
1636 Int_t AliCalorimeter::GetNvetos()
1637 {
1638 // Provide the number of veto signals associated to the calorimeter.
1639  Int_t nvetos=0;
1640  if (fVetos) nvetos=fVetos->GetEntries();
1641  return nvetos;
1642 }
1643 ///////////////////////////////////////////////////////////////////////////
1644 AliSignal* AliCalorimeter::GetVetoSignal(Int_t i)
1645 {
1646 // Provide access to the i-th veto signal of this calorimeter
1647 // Note : The first hit corresponds to i=1
1648
1649  if (i>0 && i<=GetNvetos())
1650  {
1651   return (AliSignal*)fVetos->At(i-1);
1652  }
1653  else
1654  {
1655   cout << " *AliCalorimeter::GetVetoSignal* Signal number " << i
1656        << " out of range ==> 0 returned." << endl;
1657   return 0;
1658  }
1659 }
1660 ///////////////////////////////////////////////////////////////////////////
1661 void AliCalorimeter::SetName(TString name)
1662 {
1663 // Set the name of the calorimeter system.
1664  fName=name;
1665 }
1666 ///////////////////////////////////////////////////////////////////////////
1667 TString AliCalorimeter::GetName()
1668 {
1669 // Provide the name of the calorimeter system.
1670  return fName;
1671 }
1672 ///////////////////////////////////////////////////////////////////////////
1673 void AliCalorimeter::SetSwapMode(Int_t swap)
1674 {
1675 // Set the swap mode for the module and position matrices.
1676 // At invokation of this memberfunction the default argument is swap=1.
1677 // For further details see the documentation of AliObjMatrix.
1678  if (swap==0 || swap==1)
1679  {
1680   fSwap=swap;
1681  }
1682  else
1683  {
1684   cout << " *AliCalorimeter::SetSwapMode* Invalid argument : swap = " << swap << endl; 
1685  }
1686 }
1687 ///////////////////////////////////////////////////////////////////////////
1688 Int_t AliCalorimeter::GetSwapMode()
1689 {
1690 // Provide the swap mode for the module and position matrices.
1691 // For further details see the documentation of AliObjMatrix.
1692  return fSwap;
1693 }
1694 ///////////////////////////////////////////////////////////////////////////