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