]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RALICE/AliCalorimeter.cxx
DIPO added
[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() : AliDevice()
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 }
84 ///////////////////////////////////////////////////////////////////////////
85 AliCalorimeter::~AliCalorimeter()
86 {
87 // Destructor to delete memory allocated to the various arrays and matrices
88  if (fClusters)
89  {
90   delete fClusters;
91   fClusters=0;
92  }
93  if (fVetos)
94  {
95   delete fVetos;
96   fVetos=0;
97  }
98  if (fHmodules)
99  {
100   delete fHmodules;
101   fHmodules=0;
102  }
103  if (fHclusters)
104  {
105   delete fHclusters;
106   fHclusters=0;
107  }
108  if (fMatrix)
109  {
110   delete fMatrix;
111   fMatrix=0;
112  }
113  if (fPositions)
114  {
115   delete fPositions;
116   fPositions=0;
117  }
118  if (fAttributes)
119  {
120   delete fAttributes;
121   fAttributes=0;
122  }
123 }
124 ///////////////////////////////////////////////////////////////////////////
125 AliCalorimeter::AliCalorimeter(Int_t nrow,Int_t ncol) : AliDevice()
126 {
127 // Create a calorimeter module matrix with fixed row and column size.
128 // The modules at the edges are automatically marked as "edge modules".
129  fNrows=nrow;
130  fNcolumns=ncol;
131  fClusters=0;
132
133  fSwap=0;
134  fMatrix=0;
135  fPositions=0;
136
137  fAttributes=new TObjArray(nrow);
138  fAttributes->SetOwner();
139  
140  // Mark the edge modules
141  for (Int_t row=1; row<=nrow; row++)
142  {
143   AliAttribObj* a=new AliAttribObj();
144   if (row==1 || row==nrow)
145   {
146    for (Int_t col=1; col<=ncol; col++)
147    {
148     a->SetEdgeOn(col);
149    }
150   }
151   else
152   {
153    a->SetEdgeOn(1);
154    a->SetEdgeOn(ncol);
155   }
156   fAttributes->Add(a);
157  }
158  
159  fHmodules=0;
160  fHclusters=0;
161
162  fVetos=0;
163 }
164 ///////////////////////////////////////////////////////////////////////////
165 AliCalorimeter::AliCalorimeter(const AliCalorimeter& c) : AliDevice(c)
166 {
167 // Copy constructor
168  fClusters=0;
169  fVetos=0;
170
171  fAttributes=0;
172
173  fHmodules=0;
174  fHclusters=0;
175
176  fMatrix=0;
177  fPositions=0;
178
179  fNrows=c.fNrows;
180  fNcolumns=c.fNcolumns;
181
182  fSwap=c.fSwap;
183
184  if (c.fPositions)
185  {
186   Int_t nrows=(c.fPositions)->GetMaxRow();
187   Int_t ncols=(c.fPositions)->GetMaxColumn();
188   for (Int_t irow=1; irow<=nrows; irow++)
189   {
190    for (Int_t icol=1; icol<=ncols; icol++)
191    {
192     AliPosition* p=(AliPosition*)(c.fPositions->GetObject(irow,icol));
193     if (p) SetPosition(irow,icol,*p);
194    } 
195   }  
196  }
197
198  Int_t size=0;
199  if (c.fAttributes) size=c.fAttributes->GetSize();
200  if (size)
201  {
202   fAttributes=new TObjArray(size);
203   fAttributes->SetOwner();
204   for (Int_t ia=0; ia<size; ia++)
205   {
206    AliAttribObj* a=(AliAttribObj*)(c.fAttributes->At(ia));
207    if (a) fAttributes->AddAt(new AliAttribObj(*a),ia);
208   }
209  }
210
211  Int_t n=0;
212  n=c.GetNclusters();
213  if (n)
214  {
215   fClusters=new TObjArray();
216   fClusters->SetOwner();
217   for (Int_t icl=1; icl<=n; icl++)
218   {
219    AliCalcluster* cl=c.GetCluster(icl);
220    if (cl) fClusters->Add(new AliCalcluster(*cl));
221   }
222  }
223
224  n=c.GetNvetos();
225  for (Int_t iv=1; iv<=n; iv++)
226  {
227   AliSignal* s=c.GetVetoSignal(iv);
228   if (s) AddVetoSignal(s);
229  }
230 }
231 ///////////////////////////////////////////////////////////////////////////
232 Int_t AliCalorimeter::GetNrows()
233 {
234 // Provide the number of rows for the calorimeter module matrix
235  Int_t nrows=fNrows;
236  if (!fMatrix) LoadMatrix();
237  if (fMatrix && !nrows) nrows=fMatrix->GetMaxRow();
238  return nrows;
239 }
240 ///////////////////////////////////////////////////////////////////////////
241 Int_t AliCalorimeter::GetNcolumns()
242 {
243 // Provide the number of columns for the calorimeter module matrix
244  Int_t ncols=fNcolumns;
245  if (!fMatrix) LoadMatrix();
246  if (fMatrix && !ncols) ncols=fMatrix->GetMaxColumn();
247  return ncols;
248 }
249 ///////////////////////////////////////////////////////////////////////////
250 void AliCalorimeter::SetSignal(Int_t row,Int_t col,Float_t sig)
251 {
252 // Set the signal for a certain calorimeter module.
253
254  // Check for (row,col) boundaries.
255  if (row<1 || col<1 || (fNrows && fNcolumns && (row>fNrows || col>fNcolumns)))
256  {
257   cout << " *AliCalorimeter::SetSignal* row,col : " << row << "," << col
258        << " out of range." << endl;
259   if (fNrows && fNcolumns) cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
260   return;
261  }
262
263  if (!fMatrix) LoadMatrix();
264
265  if (!fMatrix)
266  {
267   fMatrix=new AliObjMatrix();
268   fMatrix->SetSwapMode(fSwap);
269  }
270
271  AliCalmodule* mx=GetModule(row,col);
272  if (mx) // Existing module
273  {
274   mx->SetSignal(sig);
275  }
276  else // Initialise for a new module
277  {
278   AliCalmodule m;
279   m.SetRow(row);
280   m.SetColumn(col);
281   m.SetSignal(sig);
282   AliPosition* r=0;
283   if (fPositions) r=(AliPositionObj*)fPositions->GetObject(row,col);
284   if (r) m.SetPosition(*r);
285   if (fAttributes)
286   {
287    AliAttribObj* a=0;
288    if (row <= fAttributes->GetSize()) a=(AliAttribObj*)fAttributes->At(row-1);
289    if (a)
290    {
291     if (a->GetGainFlag(col)) m.SetGain(a->GetGain(col));
292     if (a->GetOffsetFlag(col)) m.SetOffset(a->GetOffset(col));
293     if (a->GetDeadValue(col)) m.SetDead();
294     if (a->GetEdgeValue(col)) m.SetEdgeValue(a->GetEdgeValue(col));
295    }
296   }
297   AddHit(m);
298   fMatrix->EnterObject(row,col,fHits->Last());
299  }
300 }
301 ///////////////////////////////////////////////////////////////////////////
302 void AliCalorimeter::AddSignal(Int_t row, Int_t col, Float_t sig)
303 {
304 // Add the signal to a certain calorimeter module.
305
306  // Check for (row,col) boundaries.
307  if (row<1 || col<1 || (fNrows && fNcolumns && (row>fNrows || col>fNcolumns)))
308  {
309   cout << " *AliCalorimeter::AddSignal* row,col : " << row << "," << col
310        << " out of range." << endl;
311   if (fNrows && fNcolumns) cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
312   return;
313  }
314  
315   AliCalmodule* m=GetModule(row,col);
316   if (!m) // initialise for new modules
317   {
318    SetSignal(row,col,sig);
319   }
320   else
321   {
322    m->AddSignal(sig);
323   }
324 }
325 ///////////////////////////////////////////////////////////////////////////
326 void AliCalorimeter::AddSignal(AliCalmodule* mod)
327 {
328 // Add the signal of module mod to the current calorimeter data.
329 // This enables mixing of calorimeter data of various events.
330 //
331 // Note : The position and attributes according to the user provided data
332 //        for the corresponding (row,col) location will be used.
333 //        In case there is no user provided data present, the position and
334 //        attributes of the first module added to the corresponding (row,col)
335 //        location will be taken, except for the "edge" and "dead" indicators.
336 //        The latter will then both be set to 0.
337
338  if (!mod) return;
339  
340  Int_t row=mod->GetRow();
341  Int_t col=mod->GetColumn();
342  Float_t sig=mod->GetSignal();
343
344  // Check for (row,col) boundaries.
345  if (row<1 || col<1 || (fNrows && fNcolumns && (row>fNrows || col>fNcolumns)))
346  {
347   cout << " *AliCalorimeter::AddSignal* row,col : " << row << "," << col
348        << " out of range." << endl;
349   if (fNrows && fNcolumns) cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
350   return;
351  }
352
353  if (!fMatrix) LoadMatrix();
354  
355  if (!fMatrix)
356  {
357   fMatrix=new AliObjMatrix();
358   fMatrix->SetSwapMode(fSwap);
359  }
360
361  AliCalmodule* mx=GetModule(row,col);
362  if (!mx) // No module existed yet at this position
363  {
364   AliCalmodule m(*mod);
365   AliPosition* r=0;
366   if (fPositions) r=(AliPositionObj*)fPositions->GetObject(row,col);
367   if (r) m.SetPosition(*r);
368   // Don't take the dead and edge attributes from this module,
369   // but from the calorimeter dbase, if present.
370   m.SetEdgeOff();
371   m.SetAlive();
372   if (fAttributes)
373   {
374    AliAttribObj* a=0;
375    if (row <= fAttributes->GetSize()) a=(AliAttribObj*)fAttributes->At(row-1);
376    if (a)
377    {
378     if (a->GetGainFlag(col)) m.SetGain(a->GetGain(col));
379     if (a->GetOffsetFlag(col)) m.SetOffset(a->GetOffset(col));
380     if (a->GetDeadValue(col)) m.SetDead();
381     if (a->GetEdgeValue(col)) m.SetEdgeValue(a->GetEdgeValue(col));
382    }
383   }
384   AddHit(m);
385   fMatrix->EnterObject(row,col,fHits->Last());
386  }
387  else
388  {
389   mx->AddSignal(sig);
390  }
391 }
392 ///////////////////////////////////////////////////////////////////////////
393 void AliCalorimeter::Reset(Int_t row,Int_t col)
394 {
395 // Reset the signal for a certain calorimeter module.
396 // Note : Module position and attributes remain unchanged.
397
398  // Check for (row,col) boundaries.
399  if (row<1 || col<1 || (fNrows && fNcolumns && (row>fNrows || col>fNcolumns)))
400  {
401   cout << " *AliCalorimeter::Reset* row,col : " << row << "," << col
402        << " out of range." << endl;
403   if (fNrows && fNcolumns) cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
404   return;
405  }
406  
407  AliCalmodule* m=GetModule(row,col);
408  if (m)
409  {
410   RemoveHit(m);
411   fMatrix->RemoveObject(row,col);
412  }
413 }
414 ///////////////////////////////////////////////////////////////////////////
415 void AliCalorimeter::Reset(Int_t mode)
416 {
417 // Reset the signals for the complete calorimeter.
418 // Normally this is done to prepare for the data of the next event.
419 //
420 // mode = 0 : Swap mode, module positions and attributes remain unchanged.
421 //        1 : Swap mode, module positions and attributes are cleared.
422 //
423 // The default is mode=0.
424 //
425 // Note : In the case of reading AliCalorimeter objects from a data file,
426 //        one has to reset the AliCalorimeter object with mode=1
427 //        (or explicitly delete it) before reading-in the next object
428 //        in order to prevent memory leaks.
429
430  if (mode<0 || mode>1)
431  {
432   cout << " *AliCalorimeter::Reset* Wrong argument. mode = " << mode << endl;
433   return;
434  }
435
436  AliDevice::Reset(mode);
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  AliPositionObj* po=0;
1144  if (fPositions) po=(AliPositionObj*)fPositions->GetObject(row,col);
1145  if (po) return po;
1146
1147  AliCalmodule* m=GetModule(row,col);
1148  return m;
1149 }
1150 ///////////////////////////////////////////////////////////////////////////
1151 Float_t AliCalorimeter::GetClusteredSignal(Int_t row,Int_t col)
1152 {
1153 // Provide the module signal after clustering.
1154
1155  // Check for (row,col) boundaries.
1156  if (row<1 || col<1 || (fNrows && fNcolumns && (row>fNrows || col>fNcolumns)))
1157  {
1158   cout << " *AliCalorimeter::GetClusteredSignal* row,col : " << row << "," << col
1159        << " out of range." << endl;
1160   if (fNrows && fNcolumns) cout << " Nrows,Ncols = " << fNrows << "," << fNcolumns << endl;
1161   return 0;
1162  }
1163  
1164  Float_t sig=0;
1165
1166  AliCalmodule* m=GetModule(row,col);
1167  if (m) sig=m->GetClusteredSignal();
1168
1169  return sig;
1170 }
1171 ///////////////////////////////////////////////////////////////////////////
1172 Int_t AliCalorimeter::GetNsignals() const
1173 {
1174 // Provide the number of modules that contain a signal
1175 // Note : The number of modules marked 'dead' but which had a signal
1176 //        are included.
1177  return GetNhits();
1178 }
1179 ///////////////////////////////////////////////////////////////////////////
1180 void AliCalorimeter::Group(Int_t n,Int_t mode)
1181 {
1182 // Group the individual modules into clusters.
1183 // Module signals of n rings around the central module will be grouped.
1184 // The grouping process will start with the module containing the highest signal
1185 // in an iterative way.
1186 // For this all fired modules are ordered w.r.t. decreasing signal.
1187 // The search mode for the module signal hierarchy can be specified by the user.
1188 //
1189 // mode = 1 : Search performed via the (row,col) structure of the matrix (SortM)
1190 //        2 : Search performed via the linear array of fired modules (SortA)
1191 //
1192 // See the docs of the memberfunctions SortM and SortA for additional details.
1193 //
1194 // Default values : n=1 mode=1.
1195
1196  if (mode<1 || mode>2)
1197  {
1198   cout << " *AliCalorimeter::Group* Invalid mode : " << mode << endl;
1199   cout << " Default value mode=1 will be used." << endl;
1200   mode=1;
1201  }
1202
1203  if (fClusters)
1204  {
1205   delete fClusters;
1206   fClusters=0;
1207  }
1208
1209  if (!fMatrix) LoadMatrix();
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   if (mode==1) SortM(); 
1220   if (mode==2) SortA();
1221
1222   Int_t nord=0;
1223   if (fOrdered) nord=fOrdered->GetEntries();
1224  
1225   // Clustering of modules. Start with the highest signal.
1226   fClusters=new TObjArray();
1227   fClusters->SetOwner();
1228   Int_t row=0;
1229   Int_t col=0;
1230   AliCalcluster* c=0;
1231   for (Int_t i=0; i<nord; i++)
1232   {
1233    AliCalmodule* m=(AliCalmodule*)fOrdered->At(i);
1234    if (!m) continue;
1235
1236    row=m->GetRow();    // row number of cluster center
1237    col=m->GetColumn(); // column number of cluster center
1238
1239    // only use modules not yet used in a cluster
1240    if (m->GetClusteredSignal() > 0.)
1241    {
1242     Int_t edge=GetEdgeValue(row,col);
1243     c=new AliCalcluster();
1244     if (!edge) c->Start(*m);   // module to start the cluster if not on edge
1245     if (c->GetNmodules() > 0)  // cluster started successfully (no edge)
1246     {
1247      fClusters->Add(c);
1248      AddRing(row,col,n); // add signals of n rings around the center
1249     }
1250     else
1251     {
1252      if (c) delete c;
1253      c=0;
1254     }
1255    }
1256   }
1257  }
1258 }
1259 ///////////////////////////////////////////////////////////////////////////
1260 void AliCalorimeter::SortM()
1261 {
1262 // Order the modules with decreasing signal by looping over the (row,col) grid
1263 // of the matrix.
1264 // Modules which were declared as "Dead" will be rejected.
1265 // The gain etc... corrected module signals will be used in the ordering process.
1266 //
1267 // Note : This method may become slow for large, very finely granulated calorimeters.
1268 //
1269 // Very specific case :
1270 // ====================
1271 // In case of various overlapping showers of which the central modules have
1272 // EXACTLY the same signal this ordering procedure may have the following
1273 // advantages and disadvantages.
1274 //
1275 // Advantages :
1276 // ------------
1277 // * In case of multi-overlapping showers, the central shower will NOT
1278 //   be "eaten-up" from both sides, resulting in a slightly more accurate
1279 //   cluster signal.
1280 // * This method produces re-producable results, irrespective of the filling
1281 //   order of the matrix modules.
1282 //
1283 // Disadvantages :
1284 // ---------------
1285 // * In case of a very high occupancy, there might be a slight effect on the
1286 //   cluster signals depending on the geometrical location in the detector matrix.
1287
1288  if (fOrdered)
1289  {
1290   delete fOrdered;
1291   fOrdered=0;
1292  }
1293
1294  Int_t nrows=fMatrix->GetMaxRow();
1295  Int_t ncols=fMatrix->GetMaxColumn();
1296
1297  Float_t signal=0.;
1298  Int_t nord=0;
1299  for (Int_t irow=1; irow<=nrows; irow++) // loop over all modules of the matrix
1300  {
1301   for (Int_t icol=1; icol<=ncols; icol++)
1302   {
1303    AliCalmodule* m=(AliCalmodule*)fMatrix->GetObject(irow,icol);
1304    if (!m) continue;
1305
1306    signal=m->GetSignal(1,1);   // get the gain etc... corrected signal
1307    if (signal <= 0.) continue; // only take alive modules with a signal
1308  
1309    if (nord == 0) // store the first module with a signal at the first ordered position
1310    {
1311     if (!fOrdered)
1312     {
1313      Int_t nhits=GetNhits();
1314      fOrdered=new TObjArray(nhits);
1315     }
1316     nord++;
1317     fOrdered->AddAt(m,nord-1);
1318     continue;
1319    }
1320  
1321    for (Int_t j=0; j<=nord; j++) // put module in the right ordered position
1322    {
1323     if (j == nord) // module has smallest signal seen so far
1324     {
1325      nord++;
1326      fOrdered->AddAt(m,j); // add module at the end
1327      break; // go for next matrix module
1328     }
1329  
1330     if (signal < ((AliCalmodule*)fOrdered->At(j))->GetSignal(1,1)) continue;
1331  
1332     nord++;
1333     for (Int_t k=nord-1; k>j; k--) // create empty position
1334     {
1335      fOrdered->AddAt(fOrdered->At(k-1),k);
1336     }
1337     fOrdered->AddAt(m,j); // put module at empty position
1338     break; // go for next matrix module
1339    }
1340   }
1341  }
1342 }
1343 ///////////////////////////////////////////////////////////////////////////
1344 void AliCalorimeter::SortA()
1345 {
1346 // Order the modules with decreasing signal by looping over the linear array
1347 // of fired modules.
1348 // Modules which were declared as "Dead" will be rejected.
1349 // The gain etc... corrected module signals will be used in the ordering process.
1350 //
1351 // Note : This method is rather fast even for large, very finely granulated calorimeters.
1352 //
1353 // Very specific case :
1354 // ====================
1355 // In case of various overlapping showers of which the central modules have
1356 // EXACTLY the same signal this ordering procedure may have the following
1357 // advantages and disadvantages.
1358 //
1359 // Advantages :
1360 // ------------
1361 // * Even in case of a very high occupancy, the resulting cluster signals
1362 //   will in general NOT depend on the geometrical location in the detector matrix.
1363 //
1364 // Disadvantages :
1365 // ---------------
1366 // * In case of multi-overlapping showers, the central shower might be
1367 //   "eaten-up" from both sides, resulting in a slightly too low value
1368 //   of the resulting cluster signal.
1369 // * This method might produce results depending on the filling
1370 //   order of the matrix modules.
1371
1372  SortHits();
1373 }
1374 ///////////////////////////////////////////////////////////////////////////
1375 void AliCalorimeter::AddRing(Int_t row, Int_t col, Int_t n)
1376 {
1377 // Add module signals of 1 ring around (row,col) to current cluster.
1378 // The gain etc... corrected module signals will be used in this process.
1379 // The parameter n denotes the maximum number of rings around cluster center.
1380 // Note : This function is used recursively.
1381
1382  if (!fMatrix) return;
1383
1384  Int_t nrows=fMatrix->GetMaxRow();
1385  Int_t ncols=fMatrix->GetMaxColumn();
1386  
1387  if (n >= 1) // Check if any rings left for recursive calls
1388  {
1389   Float_t signal=GetSignal(row,col,1); // Gain etc... corrected signal of (row,col) module
1390  
1391   Int_t lrow=row-1; if (lrow < 1) lrow=1;         // row lowerbound for ring
1392   Int_t urow=row+1; if (urow > nrows) urow=nrows; // row upperbound for ring
1393   Int_t lcol=col-1; if (lcol < 1) lcol=1;         // col lowerbound for ring
1394   Int_t ucol=col+1; if (ucol > ncols) ucol=ncols; // row upperbound for ring
1395  
1396   for (Int_t i=lrow; i<=urow; i++)
1397   {
1398    for (Int_t j=lcol; j<=ucol; j++)
1399    {
1400     // add module(i,j) to cluster if the signal <= signal(row,col)
1401     if (GetSignal(i,j,1) <= signal)
1402     {
1403      AliCalmodule* m=(AliCalmodule*)fMatrix->GetObject(i,j);
1404      if (m) ((AliCalcluster*)fClusters->At(GetNclusters()-1))->Add(*m);
1405     }
1406     AddRing(i,j,n-1); // Go for ring of modules around this (i,j) one
1407    }
1408   }
1409  }
1410 }
1411 ///////////////////////////////////////////////////////////////////////////
1412 Int_t AliCalorimeter::GetNclusters() const
1413 {
1414 // Provide the number of clusters
1415  Int_t nclu=0;
1416  if (fClusters) nclu=fClusters->GetEntries();
1417  return nclu;
1418 }
1419 ///////////////////////////////////////////////////////////////////////////
1420 AliCalcluster* AliCalorimeter::GetCluster(Int_t j) const
1421 {
1422 // Provide cluster number j
1423 // Note : j=1 denotes the first cluster
1424
1425  if (!fClusters) return 0;
1426
1427  if ((j >= 1) && (j <= GetNclusters()))
1428  {
1429   return (AliCalcluster*)fClusters->At(j-1);
1430  }
1431  else
1432  {
1433   cout << " *AliCalorimeter::GetCluster* cluster number : " << j
1434        << " out of range ==> 0 returned." << endl;
1435   return 0;
1436  }
1437 }
1438 ///////////////////////////////////////////////////////////////////////////
1439 AliCalmodule* AliCalorimeter::GetModule(Int_t j) const
1440 {
1441 // Provide 'fired' module number j
1442 // Note : j=1 denotes the first 'fired' module
1443
1444  AliCalmodule* m=(AliCalmodule*)GetHit(j);
1445  return m;
1446 }
1447 ///////////////////////////////////////////////////////////////////////////
1448 AliCalmodule* AliCalorimeter::GetModule(Int_t row,Int_t col)
1449 {
1450 // Provide access to module (row,col).
1451 // Note : first module is at (1,1).
1452
1453  AliCalmodule* m=0;
1454  if (!fMatrix) LoadMatrix();
1455  if (fMatrix) m=(AliCalmodule*)fMatrix->GetObject(row,col);
1456  return m;
1457 }
1458 ///////////////////////////////////////////////////////////////////////////
1459 TH2F* AliCalorimeter::DrawModules(Float_t thresh,Int_t mode)
1460 {
1461 // Provide a lego plot of the module signals.
1462 // The input parameter mode (default mode=0) has the same meaning as
1463 // specified in the memberfunction GetSignal(row,col,mode).
1464 // Only modules with a (corrected) signal value above the threshold
1465 // (default thresh=0) will be displayed.
1466
1467  Int_t nrows=fNrows;
1468  Int_t ncols=fNcolumns;
1469
1470  if (!fMatrix) LoadMatrix();
1471
1472  if (fMatrix && !nrows && !ncols)
1473  {
1474   nrows=fMatrix->GetMaxRow();
1475   ncols=fMatrix->GetMaxColumn();
1476  }
1477  
1478  if (fHmodules)
1479  {
1480   fHmodules->Reset();
1481  }
1482  else
1483  {
1484   fHmodules=new TH2F("fHmodules","Module signals",
1485             ncols,0.5,float(ncols)+0.5,nrows,0.5,float(nrows)+0.5);
1486  
1487   fHmodules->SetDirectory(0); // Suppress global character of histo pointer
1488  }
1489  
1490  Int_t nmods=GetNsignals();
1491
1492  Int_t row,col;
1493  Float_t signal;
1494  Int_t dead;
1495  for (Int_t i=1; i<=nmods; i++)
1496  {
1497   AliCalmodule* m=(AliCalmodule*)fMatrix->GetObject(i);
1498   if (m)
1499   {
1500    row=m->GetRow();
1501    col=m->GetColumn();
1502    dead=m->GetDeadValue();
1503    signal=0;
1504    if (!dead) signal=GetSignal(row,col,mode);
1505    if (signal>thresh) fHmodules->Fill(float(col),float(row),signal);
1506   }
1507  }
1508  
1509  fHmodules->Draw("lego");
1510  return fHmodules;
1511 }
1512 ///////////////////////////////////////////////////////////////////////////
1513 TH2F* AliCalorimeter::DrawClusters(Float_t thresh)
1514 {
1515 // Provide a lego plot of the cluster signals.
1516 // Only clusters with a signal value above the threshold (default thresh=0)
1517 // will be displayed.
1518
1519  Int_t nrows=fNrows;
1520  Int_t ncols=fNcolumns;
1521
1522  if (!fMatrix) LoadMatrix();
1523
1524  if (fMatrix && !nrows && !ncols)
1525  {
1526   nrows=fMatrix->GetMaxRow();
1527   ncols=fMatrix->GetMaxColumn();
1528  }
1529  
1530  if (fHclusters)
1531  {
1532   fHclusters->Reset();
1533  }
1534  else
1535  {
1536   fHclusters=new TH2F("fHclusters","Cluster signals",
1537             ncols,0.5,float(ncols)+0.5,nrows,0.5,float(nrows)+0.5);
1538  
1539   fHclusters->SetDirectory(0); // Suppress global character of histo pointer
1540  }
1541  
1542  AliCalcluster* c;
1543  Int_t row,col;
1544  Float_t signal;
1545  for (Int_t i=0; i<GetNclusters(); i++)
1546  {
1547   c=(AliCalcluster*)fClusters->At(i);
1548   if (c)
1549   {
1550    row=c->GetRow();
1551    col=c->GetColumn();
1552    signal=c->GetSignal();
1553    if (signal>thresh) fHclusters->Fill(float(col),float(row),signal);
1554   }
1555  }
1556  
1557  fHclusters->Draw("lego");
1558  return fHclusters;
1559 }
1560 ///////////////////////////////////////////////////////////////////////////
1561 void AliCalorimeter::Ungroup()
1562 {
1563 // Set the module signals back to the non-clustered situation
1564
1565  if (!fMatrix) LoadMatrix();
1566
1567  if (!fMatrix) return;
1568  
1569  Int_t nsig=GetNsignals();
1570
1571  Float_t signal=0;
1572  for (Int_t j=1; j<=nsig; j++)
1573  {
1574   AliCalmodule* m=(AliCalmodule*)fMatrix->GetObject(j);
1575   if (m)
1576   {
1577    signal=m->GetSignal();
1578    m->SetClusteredSignal(signal);
1579   }
1580  }
1581 }
1582 ///////////////////////////////////////////////////////////////////////////
1583 void AliCalorimeter::AddVetoSignal(AliSignal& s)
1584 {
1585 // Associate an (extrapolated) AliSignal as veto to the calorimeter.
1586  if (!fVetos)
1587  {
1588   fVetos=new TObjArray();
1589   fVetos->SetOwner();
1590  } 
1591
1592  AliSignal* sx=new AliSignal(s);
1593
1594  fVetos->Add(sx);
1595 }
1596 ///////////////////////////////////////////////////////////////////////////
1597 Int_t AliCalorimeter::GetNvetos() const
1598 {
1599 // Provide the number of veto signals associated to the calorimeter.
1600  Int_t nvetos=0;
1601  if (fVetos) nvetos=fVetos->GetEntries();
1602  return nvetos;
1603 }
1604 ///////////////////////////////////////////////////////////////////////////
1605 AliSignal* AliCalorimeter::GetVetoSignal(Int_t i) const
1606 {
1607 // Provide access to the i-th veto signal of this calorimeter
1608 // Note : The first hit corresponds to i=1
1609
1610  if (i>0 && i<=GetNvetos())
1611  {
1612   return (AliSignal*)fVetos->At(i-1);
1613  }
1614  else
1615  {
1616   cout << " *AliCalorimeter::GetVetoSignal* Signal number " << i
1617        << " out of range ==> 0 returned." << endl;
1618   return 0;
1619  }
1620 }
1621 ///////////////////////////////////////////////////////////////////////////
1622 void AliCalorimeter::SetMatrixSwapMode(Int_t swap)
1623 {
1624 // Set the swap mode for the module and position matrices.
1625 // At invokation of this memberfunction the default argument is swap=1.
1626 // For further details see the documentation of AliObjMatrix.
1627  if (swap==0 || swap==1)
1628  {
1629   fSwap=swap;
1630  }
1631  else
1632  {
1633   cout << " *AliCalorimeter::SetMatrixSwapMode* Invalid argument : swap = " << swap << endl; 
1634  }
1635 }
1636 ///////////////////////////////////////////////////////////////////////////
1637 Int_t AliCalorimeter::GetMatrixSwapMode() const
1638 {
1639 // Provide the swap mode for the module and position matrices.
1640 // For further details see the documentation of AliObjMatrix.
1641  return fSwap;
1642 }
1643 ///////////////////////////////////////////////////////////////////////////
1644 void AliCalorimeter::LoadMatrix()
1645 {
1646 // Load the matrix lookup table of module pointers from the linear hit array.
1647  Int_t nhits=GetNhits();
1648  
1649  if (!nhits) return;
1650
1651  fMatrix=new AliObjMatrix();
1652  fMatrix->SetSwapMode(fSwap);
1653
1654  Int_t row=0;
1655  Int_t col=0;
1656  for (Int_t i=1; i<=nhits; i++)
1657  {
1658   AliCalmodule* m=(AliCalmodule*)GetHit(i);
1659   if (m)
1660   {
1661    row=m->GetRow();
1662    col=m->GetColumn();
1663    fMatrix->EnterObject(row,col,m);
1664   }
1665  }
1666 }
1667 ///////////////////////////////////////////////////////////////////////////
1668 TObject* AliCalorimeter::Clone(const char* name) const
1669 {
1670 // Make a deep copy of the current object and provide the pointer to the copy.
1671 // This memberfunction enables automatic creation of new objects of the
1672 // correct type depending on the object type, a feature which may be very useful
1673 // for containers like AliEvent when adding objects in case the
1674 // container owns the objects.
1675
1676  AliCalorimeter* cal=new AliCalorimeter(*this);
1677  if (name)
1678  {
1679   if (strlen(name)) cal->SetName(name);
1680  } 
1681  return cal;
1682 }
1683 ///////////////////////////////////////////////////////////////////////////