]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONClusterFinderPeakCOG.cxx
correction of trivial typo preventing compilation
[u/mrichter/AliRoot.git] / MUON / AliMUONClusterFinderPeakCOG.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 AliMUONClusterFinderPeakCOG
20 /// 
21 /// Clusterizer class based on simple peak finder
22 ///
23 /// Pre-clustering is handled by AliMUONPreClusterFinder
24 /// From a precluster a pixel array is built, and its local maxima are used
25 /// to get pads and compute pad center of gravity.
26 ///
27 /// \author Laurent Aphecetche (for the "new" C++ structure) and 
28 /// Alexander Zinchenko, JINR Dubna, for the hardcore of it ;-)
29 //-----------------------------------------------------------------------------
30
31 #include "AliMUONClusterFinderPeakCOG.h"
32 #include "AliMUONCluster.h"
33 #include "AliMUONPad.h"
34
35 #include "AliMpPad.h"
36 #include "AliMpVSegmentation.h"
37
38 #include "AliLog.h"
39 #include "AliRunLoader.h"
40 //#include "AliCodeTimer.h"
41
42 #include <Riostream.h>
43 #include <TH2.h>
44 //#include <TCanvas.h>
45 #include <TMath.h>
46
47 /// \cond CLASSIMP
48 ClassImp(AliMUONClusterFinderPeakCOG)
49 /// \endcond
50  
51 const Double_t AliMUONClusterFinderPeakCOG::fgkZeroSuppression = 6; // average zero suppression value
52 //const Double_t AliMUONClusterFinderMLEM::fgkDistancePrecision = 1e-6; // (cm) used to check overlaps and so on
53 const Double_t AliMUONClusterFinderPeakCOG::fgkDistancePrecision = 1e-3; // (cm) used to check overlaps and so on
54 const TVector2 AliMUONClusterFinderPeakCOG::fgkIncreaseSize(-AliMUONClusterFinderPeakCOG::fgkDistancePrecision,-AliMUONClusterFinderPeakCOG::fgkDistancePrecision);
55 const TVector2 AliMUONClusterFinderPeakCOG::fgkDecreaseSize(AliMUONClusterFinderPeakCOG::fgkDistancePrecision,AliMUONClusterFinderPeakCOG::fgkDistancePrecision);
56
57 // Status flags for pads
58 const Int_t AliMUONClusterFinderPeakCOG::fgkZero = 0x0; ///< pad "basic" state
59 const Int_t AliMUONClusterFinderPeakCOG::fgkMustKeep = 0x1; ///< do not kill (for pixels)
60 const Int_t AliMUONClusterFinderPeakCOG::fgkUseForFit = 0x10; ///< should be used for fit
61 const Int_t AliMUONClusterFinderPeakCOG::fgkOver = 0x100; ///< processing is over
62 const Int_t AliMUONClusterFinderPeakCOG::fgkModified = 0x1000; ///< modified pad charge 
63 const Int_t AliMUONClusterFinderPeakCOG::fgkCoupled = 0x10000; ///< coupled pad  
64
65 //_____________________________________________________________________________
66 AliMUONClusterFinderPeakCOG::AliMUONClusterFinderPeakCOG(Bool_t plot, AliMUONVClusterFinder* clusterFinder)
67   : AliMUONVClusterFinder(),
68 fPreClusterFinder(clusterFinder),
69 fPreCluster(0x0),
70 fClusterList(),
71 fEventNumber(0),
72 fDetElemId(-1),
73 fClusterNumber(0),
74 fHistAnode(0x0),
75 fPixArray(new TObjArray(20)),
76 fDebug(0),
77 fPlot(plot),
78 fNClusters(0),
79 fNAddVirtualPads(0)
80 {
81   /// Constructor
82  
83   fSegmentation[1] = fSegmentation[0] = 0x0; 
84
85   if (fPlot) fDebug = 1;
86 }
87
88 //_____________________________________________________________________________
89 AliMUONClusterFinderPeakCOG::~AliMUONClusterFinderPeakCOG()
90 {
91 /// Destructor
92   delete fPixArray; fPixArray = 0;
93   delete fPreClusterFinder;
94   AliInfo(Form("Total clusters %d AddVirtualPad needed %d",
95                fNClusters,fNAddVirtualPads));
96 }
97
98 //_____________________________________________________________________________
99 Bool_t 
100 AliMUONClusterFinderPeakCOG::Prepare(Int_t detElemId, TClonesArray* pads[2],
101                                      const AliMpArea& area, const AliMpVSegmentation* seg[2])
102 {
103   /// Prepare for clustering
104 //  AliCodeTimerAuto("")
105   
106   for ( Int_t i = 0; i < 2; ++i )
107   {
108     fSegmentation[i] = seg[i];
109   }
110   
111   // Find out the DetElemId
112   fDetElemId = detElemId;
113   
114   // find out current event number, and reset the cluster number
115   fEventNumber = AliRunLoader::GetRunLoader()->GetEventNumber();
116   fClusterNumber = -1;
117   fClusterList.Delete();
118   
119   AliDebug(3,Form("EVT %d DE %d",fEventNumber,fDetElemId));
120   
121   if ( fPreClusterFinder->NeedSegmentation() )
122   {
123     return fPreClusterFinder->Prepare(detElemId,pads,area,seg);
124   }
125   else
126   {
127     return fPreClusterFinder->Prepare(detElemId,pads,area);
128   }
129 }
130
131 //_____________________________________________________________________________
132 AliMUONCluster* 
133 AliMUONClusterFinderPeakCOG::NextCluster()
134 {
135   /// Return next cluster
136 //  AliCodeTimerAuto("")
137   
138   // if the list of clusters is not void, pick one from there
139   TObject* o = fClusterList.At(++fClusterNumber);
140   if ( o != 0x0 ) return static_cast<AliMUONCluster*>(o);
141   
142   //FIXME : at this point, must check whether we've used all the digits
143   //from precluster : if not, let the preclustering know about those unused
144   //digits, so it can reuse them
145   
146   // if the cluster list is exhausted, we need to go to the next
147   // pre-cluster and treat it
148
149   fPreCluster = fPreClusterFinder->NextCluster();
150   
151   if (!fPreCluster)
152   {
153     // we are done
154     return 0x0;
155   }
156     
157   fClusterList.Delete(); // reset the list of clusters for this pre-cluster
158   fClusterNumber = -1; 
159   
160   WorkOnPreCluster();
161
162   // WorkOnPreCluster may have used only part of the pads, so we check that
163   // now, and let the unused pads be reused by the preclustering...
164   
165   Int_t mult = fPreCluster->Multiplicity();
166   for ( Int_t i = 0; i < mult; ++i )
167   {
168     AliMUONPad* pad = fPreCluster->Pad(i);
169     if ( !pad->IsUsed() )
170     {
171       fPreClusterFinder->UsePad(*pad);
172     }
173   }
174   
175   return NextCluster();
176 }
177
178 //_____________________________________________________________________________
179 Bool_t
180 AliMUONClusterFinderPeakCOG::WorkOnPreCluster()
181 {
182   /// Starting from a precluster, builds a pixel array, and then
183   /// extract clusters from this array
184   
185   //  AliCodeTimerAuto("")      
186
187   if (fDebug) {
188     cout << " *** Event # " << AliRunLoader::GetRunLoader()->GetEventNumber() 
189          << " det. elem.: " << fDetElemId << endl;
190     for (Int_t j = 0; j < fPreCluster->Multiplicity(); ++j) {
191       AliMUONPad* pad = fPreCluster->Pad(j);
192       printf(" bbb %3d %1d %8.4f %8.4f %8.4f %8.4f %6.1f %3d %3d %2d %1d %1d \n",
193              j, pad->Cathode(), pad->Coord(0), pad->Coord(1), pad->DX()*2, pad->DY()*2,
194              pad->Charge(), pad->Ix(), pad->Iy(), pad->Status(), pad->IsReal(), pad->IsSaturated());
195     }
196   }
197
198   AliMUONCluster* cluster = CheckPrecluster(*fPreCluster);
199   if (!cluster) return kFALSE;
200
201   BuildPixArray(*cluster);
202   
203   if ( fPixArray->GetLast() < 0 )
204   {
205     AliDebug(1,"No pixel for the above cluster");
206     delete cluster;
207     return kFALSE;
208   }
209   
210   Int_t nMax = 1, localMax[100], maxPos[100] = {0};
211   Double_t maxVal[100];
212   
213   nMax = FindLocalMaxima(fPixArray, localMax, maxVal); // find local maxima
214
215   if (nMax > 1) TMath::Sort(nMax, maxVal, maxPos, kTRUE); // in descending order
216   
217   for (Int_t i = 0; i < nMax; ++i) 
218   {
219     FindCluster(*cluster,localMax, maxPos[i]);
220   }
221
222   delete cluster;
223   if (fPlot == 0) {
224     delete fHistAnode;
225     fHistAnode = 0x0;
226   }
227   return kTRUE;
228 }
229
230 //_____________________________________________________________________________
231 Bool_t 
232 AliMUONClusterFinderPeakCOG::Overlap(const AliMUONPad& pad, const AliMUONPad& pixel)
233 {
234   /// Check if the pad and the pixel overlaps
235
236   // make a fake pad from the pixel
237   AliMUONPad tmp(pad.DetElemId(),pad.Cathode(),pad.Ix(),pad.Iy(),
238                  pixel.Coord(0),pixel.Coord(1),
239                  pixel.Size(0),pixel.Size(1),0);
240   
241   return AliMUONPad::AreOverlapping(pad,tmp,fgkDecreaseSize);
242 }
243
244 //_____________________________________________________________________________
245 AliMUONCluster* 
246 AliMUONClusterFinderPeakCOG::CheckPrecluster(const AliMUONCluster& origCluster)
247 {
248   /// Check precluster in order to attempt to simplify it (mostly for
249   /// two-cathode preclusters)
250     
251   //  AliCodeTimerAuto("")
252
253   // Disregard small clusters (leftovers from splitting or noise)
254   if ((origCluster.Multiplicity()==1 || origCluster.Multiplicity()==2) &&
255       origCluster.Charge(0)+origCluster.Charge(1) < 10) 
256   { 
257     return 0x0;
258   }
259
260   AliMUONCluster* cluster = new AliMUONCluster(origCluster);
261
262   AliDebug(2,"Start of CheckPreCluster=");
263   //StdoutToAliDebug(2,cluster->Print("full"));
264
265   AliMUONCluster* rv(0x0);
266   
267   if (cluster->Multiplicity(0) && cluster->Multiplicity(1))
268   { 
269     rv = CheckPreclusterTwoCathodes(cluster);
270   }
271   else
272   {
273     rv = cluster;
274   }
275   return rv;
276 }
277
278 //_____________________________________________________________________________
279 AliMUONCluster*
280 AliMUONClusterFinderPeakCOG::CheckPreclusterTwoCathodes(AliMUONCluster* cluster)
281 {
282   /// Check two-cathode cluster
283   
284   Int_t npad = cluster->Multiplicity();
285   Int_t* flags = new Int_t[npad];
286   for (Int_t j = 0; j < npad; ++j) flags[j] = 0;
287   
288   // Check pad overlaps
289   for ( Int_t i = 0; i < npad; ++i) 
290   {
291     AliMUONPad* padi = cluster->Pad(i);
292     if ( padi->Cathode() != 0 ) continue;
293     for (Int_t j = i+1; j < npad; ++j) 
294     {
295       AliMUONPad* padj = cluster->Pad(j);
296       if ( padj->Cathode() != 1 ) continue;
297       if ( !AliMUONPad::AreOverlapping(*padi,*padj,fgkDecreaseSize) ) continue;
298       flags[i] = flags[j] = 1; // mark overlapped pads
299     } 
300   } 
301   
302   // Check if all pads overlap
303   Int_t nFlags=0;
304   for (Int_t i = 0; i < npad; ++i) 
305   {
306     if (!flags[i]) ++nFlags;
307   }
308   
309   if (nFlags > 0) 
310   {
311     // not all pads overlap.
312     if (fDebug) cout << " nFlags: " << nFlags << endl;
313     TObjArray toBeRemoved;
314     for (Int_t i = 0; i < npad; ++i) 
315     {
316       AliMUONPad* pad = cluster->Pad(i);
317       if (flags[i]) continue;
318       Int_t cath = pad->Cathode();
319       Int_t cath1 = TMath::Even(cath);
320       // Check for edge effect (missing pads on the _other_ cathode)
321       AliMpPad mpPad = fSegmentation[cath1]->PadByPosition(pad->Position(),kFALSE);
322       if (!mpPad.IsValid()) continue;
323       //if (nFlags == 1 && pad->Charge() < fgkZeroSuppression * 3) continue;
324       if (nFlags == 1 && pad->Charge() < 20) continue;
325       AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
326                       fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),pad->Charge()));
327       toBeRemoved.AddLast(pad);
328       fPreCluster->Pad(i)->Release();
329     }
330     Int_t nRemove = toBeRemoved.GetEntriesFast();
331     for ( Int_t i = 0; i < nRemove; ++i )
332     {
333       cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
334     }
335   } 
336   
337   // Check correlations of cathode charges
338   if ( !cluster->IsSaturated() && cluster->ChargeAsymmetry() > 1 )
339   {
340     // big difference
341     Int_t cathode = cluster->MaxRawChargeCathode();
342     Int_t imin(-1);
343     Int_t imax(-1);
344     Double_t cmax(0);
345     Double_t cmin(1E9);
346     
347     // get min and max pad charges on the cathode opposite to the 
348     // max pad (given by MaxRawChargeCathode())
349     //
350     Int_t mult = cluster->Multiplicity();
351     for ( Int_t i = 0; i < mult; ++i )
352     {
353       AliMUONPad* pad = cluster->Pad(i);
354       if ( pad->Cathode() != cathode || !pad->IsReal() )
355       {
356         // only consider pads in the opposite cathode, and
357         // only consider real pads (i.e. exclude the virtual ones)
358         continue;
359       }
360       if ( pad->Charge() < cmin )
361       {
362         cmin = pad->Charge();
363         imin = i;
364         if (imax < 0) {
365           imax = imin;
366           cmax = cmin;
367         }
368       }
369       else if ( pad->Charge() > cmax )
370       {
371         cmax = pad->Charge();
372         imax = i;
373       }      
374     }
375     AliDebug(2,Form("Pad imin,imax %d,%d cmin,cmax %e,%e",
376                     imin,imax,cmin,cmax));
377     //
378     // arrange pads according to their distance to the max, normalized
379     // to the pad size
380     Double_t* dist = new Double_t[mult];
381     Double_t dxMin(1E9);
382     Double_t dyMin(1E9);
383     Double_t dmin(0);
384     
385     AliMUONPad* padmax = cluster->Pad(imax);
386     
387     for ( Int_t i = 0; i < mult; ++i )
388     {
389       dist[i] = 0.0;
390       if ( i == imax) continue;
391       AliMUONPad* pad = cluster->Pad(i);
392       if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
393       Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
394       Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
395       dist[i] = TMath::Sqrt(dx*dx+dy*dy);
396       if ( i == imin )
397       {
398         dmin = dist[i] + 1E-3; // distance to the pad with minimum charge
399         dxMin = dx;
400         dyMin = dy;
401       }      
402     }
403     
404     TMath::Sort(mult,dist,flags,kFALSE); // in ascending order
405     Double_t xmax(-1), distPrev(999);
406     TObjArray toBeRemoved;
407     
408     for ( Int_t i = 0; i < mult; ++i )
409     {
410       Int_t indx = flags[i];
411       AliMUONPad* pad = cluster->Pad(indx);
412       if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
413       if ( dist[indx] > dmin )
414       {
415         // farther than the minimum pad
416         Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
417         Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
418         dx *= dxMin;
419         dy *= dyMin;
420         if (dx >= 0 && dy >= 0) continue;
421         if (TMath::Abs(dx) > TMath::Abs(dy) && dx >= 0) continue;
422         if (TMath::Abs(dy) > TMath::Abs(dx) && dy >= 0) continue;        
423       }
424       if (dist[indx] > distPrev + 1) break; // overstepping empty pads
425       if ( pad->Charge() <= cmax || TMath::Abs(dist[indx]-xmax) < 1E-3 )
426       {
427         // release pad
428         if (TMath::Abs(dist[indx]-xmax) < 1.e-3) 
429         {
430           cmax = TMath::Max(pad->Charge(),cmax);
431         }
432         else
433         {
434           cmax = pad->Charge();
435         }
436         xmax = dist[indx];
437         distPrev = dist[indx];
438         AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
439                         fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),
440                         pad->Charge()));
441   
442         toBeRemoved.AddLast(pad);
443         fPreCluster->Pad(indx)->Release();
444       }
445     }
446     Int_t nRemove = toBeRemoved.GetEntriesFast();
447     for ( Int_t i = 0; i < nRemove; ++i )
448     {
449       cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
450     }    
451     delete[] dist;
452   } // if ( !cluster->IsSaturated() && 
453   
454   delete[] flags;
455   
456   AliDebug(2,"End of CheckPreClusterTwoCathodes=");
457   //StdoutToAliDebug(2,cluster->Print("full"));
458
459   return cluster;    
460 }
461
462 //_____________________________________________________________________________
463 void
464 AliMUONClusterFinderPeakCOG::CheckOverlaps()
465 {
466   /// For debug only : check if some pixels overlap...
467   
468   Int_t nPix = fPixArray->GetLast()+1;
469   Int_t dummy(0);
470   
471   for ( Int_t i = 0; i < nPix; ++i )
472   {
473     AliMUONPad* pixelI = Pixel(i);
474     AliMUONPad pi(dummy,dummy,dummy,dummy,
475                   pixelI->Coord(0),pixelI->Coord(1),
476                   pixelI->Size(0),pixelI->Size(1),0.0);
477     
478     for ( Int_t j = i+1; j < nPix; ++j )
479     {
480       AliMUONPad* pixelJ = Pixel(j);
481       AliMUONPad pj(dummy,dummy,dummy,dummy,
482                     pixelJ->Coord(0),pixelJ->Coord(1),
483                     pixelJ->Size(0),pixelJ->Size(1),0.0);  
484       AliMpArea area;
485       
486       if ( AliMUONPad::AreOverlapping(pi,pj,fgkDecreaseSize,area) )
487       {
488         AliInfo(Form("The following 2 pixels (%d and %d) overlap !",i,j));
489         /*
490         StdoutToAliInfo(pixelI->Print();
491                         cout << " Surface = " << pixelI->Size(0)*pixelI->Size(1)*4 << endl;
492                         pixelJ->Print();
493                         cout << " Surface = " << pixelJ->Size(0)*pixelJ->Size(1)*4 << endl;
494                         cout << " Area surface = " << area.Dimensions().X()*area.Dimensions().Y()*4 << endl;
495                         cout << "-------" << endl;
496                         );
497         */        
498       }
499     }    
500   }
501 }
502
503 //_____________________________________________________________________________
504 void AliMUONClusterFinderPeakCOG::BuildPixArray(AliMUONCluster& cluster)
505 {
506   /// Build pixel array 
507   
508   Int_t npad = cluster.Multiplicity();
509   if (npad<=0) 
510   {
511     AliWarning("Got no pad at all ?!");
512   }
513   
514   fPixArray->Delete();
515   BuildPixArrayOneCathode(cluster);
516   
517 //  StdoutToAliDebug(2,cout << "End of BuildPixelArray:" << endl;
518 //                   fPixArray->Print(););
519   //CheckOverlaps();//FIXME : this is for debug only. Remove it.
520 }
521
522 //_____________________________________________________________________________
523 void AliMUONClusterFinderPeakCOG::BuildPixArrayOneCathode(AliMUONCluster& cluster)
524 {
525   /// Build the pixel array
526
527 //  AliDebug(2,Form("cluster.Multiplicity=%d",cluster.Multiplicity()));
528
529   TVector2 dim = cluster.MinPadDimensions (-1, kFALSE);
530   Double_t width[2] = {dim.X(), dim.Y()}, xy0[2];
531   Int_t found[2] = {0}, mult = cluster.Multiplicity();
532
533   for ( Int_t i = 0; i < mult; ++i) {
534     AliMUONPad* pad = cluster.Pad(i);
535     for (Int_t j = 0; j < 2; ++j) {
536       if (found[j] == 0 && TMath::Abs(pad->Size(j)-width[j]) < fgkDistancePrecision) { 
537         xy0[j] = pad->Coord(j);
538         found[j] = 1;
539       }
540     }
541     if (found[0] && found[1]) break;
542   }
543
544   Double_t min[2], max[2];
545   Int_t cath0 = 0, cath1 = 1;
546   if (cluster.Multiplicity(0) == 0) cath0 = 1;
547   else if (cluster.Multiplicity(1) == 0) cath1 = 0;
548
549   TVector2 leftDown = cluster.Area(cath0).LeftDownCorner();
550   TVector2 rightUp = cluster.Area(cath0).RightUpCorner();
551   min[0] = leftDown.X();
552   min[1] = leftDown.Y();
553   max[0] = rightUp.X();
554   max[1] = rightUp.Y();
555   if (cath1 != cath0) {
556     leftDown = cluster.Area(cath1).LeftDownCorner();
557     rightUp = cluster.Area(cath1).RightUpCorner();
558     min[0] = TMath::Max (min[0], leftDown.X());
559     min[1] = TMath::Max (min[1], leftDown.Y());
560     max[0] = TMath::Min (max[0], rightUp.X());
561     max[1] = TMath::Min (max[1], rightUp.Y());
562   }
563
564   // Adjust limits
565   //width[0] /= 2; width[1] /= 2; // just for check
566   Int_t nbins[2];
567   for (Int_t i = 0; i < 2; ++i) {
568     Double_t dist = (min[i] - xy0[i]) / width[i] / 2;
569     if (TMath::Abs(dist) < 1.e-6) dist = -1.e-6;
570     min[i] = xy0[i] + (TMath::Nint(dist-TMath::Sign(1.e-6,dist)) 
571                        + TMath::Sign(0.5,dist)) * width[i] * 2;
572     nbins[i] = TMath::Nint ((max[i] - min[i]) / width[i] / 2);
573     if (nbins[i] == 0) ++nbins[i];
574     max[i] = min[i] + nbins[i] * width[i] * 2;
575     //cout << dist << " " << min[i] << " " << max[i] << " " << nbins[i] << endl;
576   }
577
578   // Book histogram
579   TH2D *hist1 = new TH2D ("Grid", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
580   TH2D *hist2 = new TH2D ("Entries", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
581   TAxis *xaxis = hist1->GetXaxis();
582   TAxis *yaxis = hist1->GetYaxis();
583
584   // Fill histogram
585   for ( Int_t i = 0; i < mult; ++i) {
586     AliMUONPad* pad = cluster.Pad(i);
587     Int_t ix0 = xaxis->FindBin(pad->X());
588     Int_t iy0 = yaxis->FindBin(pad->Y());
589     PadOverHist(0, ix0, iy0, pad, hist1, hist2);
590   }
591
592   // Store pixels
593   for (Int_t i = 1; i <= nbins[0]; ++i) {
594     Double_t x = xaxis->GetBinCenter(i);
595     for (Int_t j = 1; j <= nbins[1]; ++j) {
596       if (hist2->GetCellContent(i,j) < 0.1) continue;
597       if (cath0 != cath1) {
598         // Two-sided cluster
599         Double_t cont = hist2->GetCellContent(i,j);
600         if (cont < 999.) continue;
601         if (cont-Int_t(cont/1000.)*1000. < 0.5) continue;
602       }
603       //if (hist2->GetCellContent(i,j) < 1.1 && cluster.Multiplicity(0) && 
604       //  cluster.Multiplicity(1)) continue;
605       Double_t y = yaxis->GetBinCenter(j);
606       Double_t charge = hist1->GetCellContent(i,j);
607       AliMUONPad* pixPtr = new AliMUONPad(x, y, width[0], width[1], charge);
608       fPixArray->Add(pixPtr);
609     }  
610   }
611   /*
612   if (fPixArray->GetEntriesFast() == 1) {
613     // Split pixel into 2
614     AliMUONPad* pixPtr = static_cast<AliMUONPad*> (fPixArray->UncheckedAt(0));
615     pixPtr->SetSize(0,width[0]/2.);
616     pixPtr->Shift(0,-width[0]/4.);
617     pixPtr = new AliMUONPad(pixPtr->X()+width[0], pixPtr->Y(), width[0]/2., width[1], pixPtr->Charge());
618     fPixArray->Add(pixPtr);
619   }
620   */
621   //fPixArray->Print();
622   delete hist1;
623   delete hist2;
624 }
625
626 //_____________________________________________________________________________
627 void AliMUONClusterFinderPeakCOG::PadOverHist(Int_t idir, Int_t ix0, Int_t iy0, AliMUONPad *pad,
628                                               TH2D *hist1, TH2D *hist2)
629 {
630   /// "Span" pad over histogram in the direction idir
631
632   TAxis *axis = idir == 0 ? hist1->GetXaxis() : hist1->GetYaxis();
633   Int_t nbins = axis->GetNbins(), cath = pad->Cathode();
634   Double_t bin = axis->GetBinWidth(1), amask = TMath::Power(1000.,cath*1.);
635
636   Int_t nbinPad = (Int_t)(pad->Size(idir)/bin*2+fgkDistancePrecision) + 1; // number of bins covered by pad
637
638   for (Int_t i = 0; i < nbinPad; ++i) {
639     Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
640     if (ixy > nbins) break;
641     Double_t lowEdge = axis->GetBinLowEdge(ixy);
642     if (lowEdge + fgkDistancePrecision > pad->Coord(idir) + pad->Size(idir)) break;
643     if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
644     else {
645       // Fill histogram
646       Double_t cont = pad->Charge();
647       if (hist2->GetCellContent(ix0, ixy) > 0.1) 
648         cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont) 
649              + TMath::Min (TMath::Max(hist1->GetCellContent(ix0, ixy),cont)*0.1, 10.);
650       hist1->SetCellContent(ix0, ixy, cont);
651       hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
652     }
653   }
654
655   for (Int_t i = -1; i > -nbinPad; --i) {
656     Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
657     if (ixy < 1) break;
658     Double_t upEdge = axis->GetBinUpEdge(ixy);
659     if (upEdge - fgkDistancePrecision < pad->Coord(idir) - pad->Size(idir)) break;
660     if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
661     else {
662       // Fill histogram
663       Double_t cont = pad->Charge();
664       if (hist2->GetCellContent(ix0, ixy) > 0.1) 
665         cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont)
666              + TMath::Min (TMath::Max(hist1->GetCellContent(ix0, ixy),cont)*0.1,10.);
667       hist1->SetCellContent(ix0, ixy, cont);
668       hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
669     }
670   }
671 }
672
673 //_____________________________________________________________________________
674 Int_t AliMUONClusterFinderPeakCOG::FindLocalMaxima(TObjArray *pixArray, Int_t *localMax, Double_t *maxVal)
675 {
676 /// Find local maxima in pixel space 
677
678   AliDebug(1,Form("nPix=%d",pixArray->GetLast()+1));
679
680   //TH2D *hist = NULL;
681   //delete ((TH2D*) gROOT->FindObject("anode"));
682   //if (pixArray == fPixArray) hist = (TH2D*) gROOT->FindObject("anode");
683   //else { hist = (TH2D*) gROOT->FindObject("anode1"); cout << hist << endl; }
684   //if (hist) hist->Delete();
685   delete fHistAnode;
686  
687   Double_t xylim[4] = {999, 999, 999, 999};
688
689   Int_t nPix = pixArray->GetEntriesFast();
690   AliMUONPad *pixPtr = 0;
691   for (Int_t ipix = 0; ipix < nPix; ++ipix) {
692     pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
693     for (Int_t i = 0; i < 4; ++i) 
694          xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
695   }
696   for (Int_t i = 0; i < 4; ++i) xylim[i] -= pixPtr->Size(i/2); 
697
698   Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
699   Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
700   if (pixArray == fPixArray) fHistAnode = new TH2D("anode","anode",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
701   else fHistAnode = new TH2D("anode1","anode1",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
702   for (Int_t ipix = 0; ipix < nPix; ++ipix) {
703     pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
704     fHistAnode->Fill(pixPtr->Coord(0), pixPtr->Coord(1), pixPtr->Charge());
705   }
706 //  if (fDraw && pixArray == fPixArray) fDraw->DrawHist("c2", hist);
707
708   Int_t nMax = 0, indx, nxy = ny * nx;
709   Int_t *isLocalMax = new Int_t[nxy];
710   for (Int_t i = 0; i < nxy; ++i) isLocalMax[i] = 0; 
711
712   for (Int_t i = 1; i <= ny; ++i) {
713     indx = (i-1) * nx;
714     for (Int_t j = 1; j <= nx; ++j) {
715       if (fHistAnode->GetCellContent(j,i) < 0.5) continue;
716       //if (isLocalMax[indx+j-1] < 0) continue;
717       if (isLocalMax[indx+j-1] != 0) continue;
718       FlagLocalMax(fHistAnode, i, j, isLocalMax);
719     }
720   }
721
722   for (Int_t i = 1; i <= ny; ++i) {
723     indx = (i-1) * nx;
724     for (Int_t j = 1; j <= nx; ++j) {
725       if (isLocalMax[indx+j-1] > 0) { 
726         localMax[nMax] = indx + j - 1; 
727         maxVal[nMax++] = fHistAnode->GetCellContent(j,i);
728         if (nMax > 99) break;
729       }
730     }
731     if (nMax > 99) {
732       AliError(" Too many local maxima !!!");
733       break;
734     }
735   }
736   if (fDebug) cout << " Local max: " << nMax << endl;
737   delete [] isLocalMax; 
738   return nMax;
739 }
740
741 //_____________________________________________________________________________
742 void AliMUONClusterFinderPeakCOG::FlagLocalMax(TH2D *hist, Int_t i, Int_t j, Int_t *isLocalMax)
743 {
744 /// Flag pixels (whether or not local maxima)
745
746   Int_t nx = hist->GetNbinsX();
747   Int_t ny = hist->GetNbinsY();
748   Int_t cont = TMath::Nint (hist->GetCellContent(j,i));
749   Int_t cont1 = 0, indx = (i-1)*nx+j-1, indx1 = 0, indx2 = 0;
750
751   Int_t ie = i + 2, je = j + 2;
752   for (Int_t i1 = i-1; i1 < ie; ++i1) {
753     if (i1 < 1 || i1 > ny) continue;
754     indx1 = (i1 - 1) * nx;
755     for (Int_t j1 = j-1; j1 < je; ++j1) {
756       if (j1 < 1 || j1 > nx) continue;
757       if (i == i1 && j == j1) continue;
758       indx2 = indx1 + j1 - 1;
759       cont1 = TMath::Nint (hist->GetCellContent(j1,i1));
760       if (cont < cont1) { isLocalMax[indx] = -1; return; }
761       else if (cont > cont1) isLocalMax[indx2] = -1;
762       else { // the same charge
763         isLocalMax[indx] = 1; 
764         if (isLocalMax[indx2] == 0) {
765           FlagLocalMax(hist, i1, j1, isLocalMax);
766           if (isLocalMax[indx2] < 0) { isLocalMax[indx] = -1; return; }
767           else isLocalMax[indx2] = -1;
768         }
769       } 
770     }
771   }
772   isLocalMax[indx] = 1; // local maximum
773 }
774
775 //_____________________________________________________________________________
776 void AliMUONClusterFinderPeakCOG::FindCluster(AliMUONCluster& cluster, 
777                                            Int_t *localMax, Int_t iMax)
778 {
779 /// Find pixel cluster around local maximum \a iMax and pick up pads
780 /// overlapping with it
781
782   //TH2D *hist = (TH2D*) gROOT->FindObject("anode");
783   /* Just for check
784   TCanvas* c = new TCanvas("Anode","Anode",800,600);
785   c->cd();
786   hist->Draw("lego1Fb"); // debug
787   c->Update();
788   Int_t tmp;
789   cin >> tmp;
790   */
791   Int_t nx = fHistAnode->GetNbinsX();
792   //Int_t ny = hist->GetNbinsY();
793   Int_t ic = localMax[iMax] / nx + 1;
794   Int_t jc = localMax[iMax] % nx + 1;
795
796   // Get min pad dimensions for the precluster
797   Int_t nSides = 2;
798   if (cluster.Multiplicity(0) == 0 || cluster.Multiplicity(1) == 0) nSides = 1;
799   TVector2 dim0 = cluster.MinPadDimensions(0, -1, kFALSE);
800   TVector2 dim1 = cluster.MinPadDimensions(1, -1, kFALSE);
801   //Double_t width[2][2] = {{dim0.X(), dim0.Y()},{dim1.X(),dim1.Y()}};
802   Int_t nonb[2] = {1, 0}; // coordinate index vs cathode
803   if (nSides == 1 || dim0.X() < dim1.X() - fgkDistancePrecision) {
804     nonb[0] = 0;
805     nonb[1] = 1;
806   }
807   Bool_t samex = kFALSE, samey = kFALSE;
808   if (TMath::Abs(dim0.X()-dim1.X()) < fgkDistancePrecision) samex = kTRUE; // the same X pad size on both planes 
809   if (TMath::Abs(dim0.Y()-dim1.Y()) < fgkDistancePrecision) samey = kTRUE; // the same Y pad size on both planes 
810
811   // Drop all pixels from the array - pick up only the ones from the cluster
812   //fPixArray->Delete();
813
814   Double_t wx = fHistAnode->GetXaxis()->GetBinWidth(1)/2; 
815   Double_t wy = fHistAnode->GetYaxis()->GetBinWidth(1)/2;  
816   Double_t yc = fHistAnode->GetYaxis()->GetBinCenter(ic);
817   Double_t xc = fHistAnode->GetXaxis()->GetBinCenter(jc);
818   Double_t cont = fHistAnode->GetCellContent(jc,ic);
819   AliMUONPad* pixPtr = new AliMUONPad (xc, yc, wx, wy, cont);
820   if (fDebug) pixPtr->Print("full"); 
821
822   Int_t npad = cluster.Multiplicity();
823   
824   // Pick up pads which overlap with the maximum pixel and find pads with the max signal
825   Double_t qMax[2] = {0}; 
826   AliMUONPad *matrix[2][3] = {{0x0,0x0,0x0},{0x0,0x0,0x0}};
827   for (Int_t j = 0; j < npad; ++j) 
828   {
829     AliMUONPad* pad = cluster.Pad(j);
830     if ( Overlap(*pad,*pixPtr) )
831     {
832       if (fDebug) { cout << j << " "; pad->Print("full"); }
833       if (pad->Charge() > qMax[pad->Cathode()]) {
834         qMax[pad->Cathode()] = pad->Charge();
835         matrix[pad->Cathode()][1] = pad;
836         if (nSides == 1) matrix[!pad->Cathode()][1] = pad;
837       }
838     }
839   }
840   //if (nSides == 2 && (matrix[0][1] == 0x0 || matrix[1][1] == 0x0)) return; // ???
841
842   // Find neighbours of maxima to have 3 pads per direction (if possible)
843   for (Int_t j = 0; j < npad; ++j) 
844   {
845     AliMUONPad* pad = cluster.Pad(j);
846     Int_t cath = pad->Cathode();
847     if (pad == matrix[cath][1]) continue;
848     Int_t nLoops = 3 - nSides;
849
850     for (Int_t k = 0; k < nLoops; ++k) {
851       Int_t cath1 = cath;
852       if (k) cath1 = !cath;
853
854       // Check the coordinate corresponding to the cathode (bending or non-bending case)
855       Double_t dist = pad->Coord(nonb[cath1]) - matrix[cath][1]->Coord(nonb[cath1]);
856       Double_t dir = TMath::Sign (1., dist);
857       dist = TMath::Abs(dist) - pad->Size(nonb[cath1]) - matrix[cath][1]->Size(nonb[cath1]);
858
859       if (TMath::Abs(dist) < fgkDistancePrecision) {
860         // Check the other coordinate
861         dist = pad->Coord(!nonb[cath1]) - matrix[cath1][1]->Coord(!nonb[cath1]);
862         if (TMath::Abs(dist) > 
863             TMath::Max(pad->Size(!nonb[cath1]), matrix[cath1][1]->Size(!nonb[cath1])) - fgkDistancePrecision) break;
864         Int_t idir = TMath::Nint (dir);
865         if (matrix[cath1][1+idir] == 0x0) matrix[cath1][1+idir] = pad;
866         else if (pad->Charge() > matrix[cath1][1+idir]->Charge()) matrix[cath1][1+idir] = pad; // diff. segmentation
867         //cout << pad->Coord(nonb[cath1]) << " " << pad->Coord(!nonb[cath1]) << " " << pad->Size(nonb[cath1]) << " " << pad->Size(!nonb[cath1]) << " " << pad->Charge() << endl ;
868         break;
869       }
870     }
871   }
872
873   Double_t coord[2] = {0.}, qAver = 0.;
874   for (Int_t i = 0; i < 2; ++i) {
875     Double_t q = 0.;
876     Double_t coordQ = 0.;
877     Int_t cath = matrix[i][1]->Cathode();
878     if (i && nSides == 1) cath = !cath;
879     for (Int_t j = 0; j < 3; ++j) {
880       if (matrix[i][j] == 0x0) continue;
881       Double_t dq = matrix[i][j]->Charge();
882       q += dq;
883       coordQ += dq * matrix[i][j]->Coord(nonb[cath]);
884       //coordQ += (matrix[i][j]->Charge() * matrix[i][j]->Coord(nonb[cath]));
885     }
886     coord[cath] = coordQ / q;
887     qAver = TMath::Max (qAver, q);
888   }
889
890   //qAver = TMath::Sqrt(qAver);
891   if ( qAver >= 14 )
892   {
893     
894     AliMUONCluster* cluster1 = new AliMUONCluster(cluster);
895       
896     cluster1->SetCharge(qAver,qAver);
897     if (nonb[0] == 1) 
898       cluster1->SetPosition(TVector2(coord[1],coord[0]),TVector2(0.,0.));
899     else 
900       cluster1->SetPosition(TVector2(coord[0],coord[1]),TVector2(0.,0.));
901
902     cluster1->SetChi2(0.);
903       
904     // FIXME: we miss some information in this cluster, as compared to 
905     // the original AddRawCluster code.
906       
907     AliDebug(2,Form("Adding RawCluster detElemId %4d mult %2d charge %5d (xl,yl)=(%9.6g,%9.6g)",
908                     fDetElemId,cluster1->Multiplicity(),(Int_t)cluster1->Charge(),
909                     cluster1->Position().X(),cluster1->Position().Y()));
910         
911     fClusterList.Add(cluster1);
912   }
913 }
914
915 //_____________________________________________________________________________
916 AliMUONClusterFinderPeakCOG&  
917 AliMUONClusterFinderPeakCOG::operator=(const AliMUONClusterFinderPeakCOG& rhs)
918 {
919 /// Protected assignement operator
920
921   if (this == &rhs) return *this;
922
923   AliFatal("Not implemented.");
924     
925   return *this;  
926 }    
927
928 //_____________________________________________________________________________
929 void AliMUONClusterFinderPeakCOG::PadsInXandY(AliMUONCluster& cluster,
930                                            Int_t &nInX, Int_t &nInY) const
931 {
932   /// Find number of pads in X and Y-directions (excluding virtual ones and
933   /// overflows)
934
935   //Int_t statusToTest = 1;
936   Int_t statusToTest = fgkUseForFit;
937   
938   //if ( nInX < 0 ) statusToTest = 0;
939   if ( nInX < 0 ) statusToTest = fgkZero;
940        
941   Bool_t mustMatch(kTRUE);
942
943   AliMpIntPair cn = cluster.NofPads(statusToTest,mustMatch);
944   
945   nInX = cn.GetFirst();
946   nInY = cn.GetSecond();
947 }
948
949 //_____________________________________________________________________________
950 void AliMUONClusterFinderPeakCOG::RemovePixel(Int_t i)
951 {
952   /// Remove pixel at index i
953   AliMUONPad* pixPtr = Pixel(i);
954   fPixArray->RemoveAt(i); 
955   delete pixPtr;
956 }
957
958 //_____________________________________________________________________________
959 AliMUONPad* 
960 AliMUONClusterFinderPeakCOG::Pixel(Int_t i) const
961 {
962   /// Returns pixel at index i
963   return static_cast<AliMUONPad*>(fPixArray->UncheckedAt(i));
964 }
965
966 //_____________________________________________________________________________
967 void 
968 AliMUONClusterFinderPeakCOG::Print(Option_t* what) const
969 {
970   /// printout
971   TString swhat(what);
972   swhat.ToLower();
973   if ( swhat.Contains("precluster") )
974   {
975     if ( fPreCluster) fPreCluster->Print();
976   }
977 }
978
979