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