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