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