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