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