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