]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONClusterFinderMLEM.cxx
In MUONRecoCheck:
[u/mrichter/AliRoot.git] / MUON / AliMUONClusterFinderMLEM.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 AliMUONClusterFinderMLEM
20 /// 
21 /// Clusterizer class based on the Expectation-Maximization algorithm
22 ///
23 /// Pre-clustering is handled by AliMUONPreClusterFinder
24 /// From a precluster a pixel array is built, and from this pixel array
25 /// a list of clusters is output (using AliMUONClusterSplitterMLEM).
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 "AliMUONClusterFinderMLEM.h"
32 #include "AliLog.h"
33 #include "AliMUONCluster.h"
34 #include "AliMUONClusterSplitterMLEM.h"
35 #include "AliMUONVDigit.h"
36 #include "AliMUONPad.h"
37 #include "AliMUONPreClusterFinder.h"
38 #include "AliMpPad.h"
39 #include "AliMpVPadIterator.h"
40 #include "AliMpVSegmentation.h"
41 #include "AliRunLoader.h"
42 #include "AliMUONVDigitStore.h"
43 #include <Riostream.h>
44 #include <TH2.h>
45 #include <TMinuit.h>
46 #include <TCanvas.h>
47 #include <TMath.h>
48 #include "AliCodeTimer.h"
49
50 /// \cond CLASSIMP
51 ClassImp(AliMUONClusterFinderMLEM)
52 /// \endcond
53  
54 const Double_t AliMUONClusterFinderMLEM::fgkDistancePrecision = 1e-3; // (cm) used to check overlaps and so on
55 const TVector2 AliMUONClusterFinderMLEM::fgkIncreaseSize(-AliMUONClusterFinderMLEM::fgkDistancePrecision,-AliMUONClusterFinderMLEM::fgkDistancePrecision);
56 const TVector2 AliMUONClusterFinderMLEM::fgkDecreaseSize(AliMUONClusterFinderMLEM::fgkDistancePrecision,AliMUONClusterFinderMLEM::fgkDistancePrecision);
57
58 // Status flags for pads
59 const Int_t AliMUONClusterFinderMLEM::fgkZero = 0x0; ///< pad "basic" state
60 const Int_t AliMUONClusterFinderMLEM::fgkMustKeep = 0x1; ///< do not kill (for pixels)
61 const Int_t AliMUONClusterFinderMLEM::fgkUseForFit = 0x10; ///< should be used for fit
62 const Int_t AliMUONClusterFinderMLEM::fgkOver = 0x100; ///< processing is over
63 const Int_t AliMUONClusterFinderMLEM::fgkModified = 0x1000; ///< modified pad charge 
64 const Int_t AliMUONClusterFinderMLEM::fgkCoupled = 0x10000; ///< coupled pad  
65
66 //_____________________________________________________________________________
67 AliMUONClusterFinderMLEM::AliMUONClusterFinderMLEM(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 fHistMlem(0x0),
76 fHistAnode(0x0),
77 fPixArray(new TObjArray(20)),
78 fDebug(0),
79 fPlot(plot),
80 fSplitter(0x0),
81 fNClusters(0),
82 fNAddVirtualPads(0),
83 fLowestPixelCharge(0),
84 fLowestPadCharge(0),
85 fLowestClusterCharge(0)
86 {
87   /// Constructor
88
89   fkSegmentation[1] = fkSegmentation[0] = 0x0; 
90
91   if (fPlot) fDebug = 1;
92 }
93
94 //_____________________________________________________________________________
95 AliMUONClusterFinderMLEM::~AliMUONClusterFinderMLEM()
96 {
97 /// Destructor
98   delete fPixArray; fPixArray = 0;
99 //  delete fDraw;
100   delete fPreClusterFinder;
101   delete fSplitter;
102   AliInfo(Form("Total clusters %d AddVirtualPad needed %d",
103                fNClusters,fNAddVirtualPads));
104 }
105
106 //_____________________________________________________________________________
107 Bool_t 
108 AliMUONClusterFinderMLEM::Prepare(Int_t detElemId,
109                                   TClonesArray* pads[2],
110                                   const AliMpArea& area,
111                                   const AliMpVSegmentation* seg[2])
112 {
113   /// Prepare for clustering
114 //  AliCodeTimerAuto("",0)
115   
116   for ( Int_t i = 0; i < 2; ++i )
117   {
118     fkSegmentation[i] = seg[i];
119   }
120   
121   // Find out the DetElemId
122   fDetElemId = detElemId;
123   
124   delete fSplitter;
125   fSplitter = new AliMUONClusterSplitterMLEM(fDetElemId,
126                                              fPixArray,
127                                              fLowestPixelCharge,
128                                              fLowestPadCharge,
129                                              fLowestClusterCharge);
130   fSplitter->SetDebug(fDebug);
131     
132   // find out current event number, and reset the cluster number
133   AliRunLoader *runLoader = AliRunLoader::Instance();
134   fEventNumber = runLoader ? runLoader->GetEventNumber() : 0;
135   fClusterNumber = -1;
136   fClusterList.Delete();
137   fPixArray->Delete();
138
139   AliDebug(3,Form("EVT %d DE %d",fEventNumber,fDetElemId));
140   
141   if ( fPreClusterFinder->NeedSegmentation() )
142   {
143     return fPreClusterFinder->Prepare(detElemId,pads,area,seg);
144   }
145   else
146   {
147     return fPreClusterFinder->Prepare(detElemId,pads,area);
148   }
149 }
150
151 //_____________________________________________________________________________
152 AliMUONCluster* 
153 AliMUONClusterFinderMLEM::NextCluster()
154 {
155   /// Return next cluster
156 //  AliCodeTimerAuto("",0)
157   
158   // if the list of clusters is not void, pick one from there
159   TObject* o = fClusterList.At(++fClusterNumber);
160   if ( o != 0x0 ) return static_cast<AliMUONCluster*>(o);
161   
162   //FIXME : at this point, must check whether we've used all the digits
163   //from precluster : if not, let the preclustering know about those unused
164   //digits, so it can reuse them
165   
166   // if the cluster list is exhausted, we need to go to the next
167   // pre-cluster and treat it
168
169   fPreCluster = fPreClusterFinder->NextCluster();
170
171   fPixArray->Delete();
172   fClusterList.Delete(); // reset the list of clusters for this pre-cluster
173   fClusterNumber = -1; //AZ
174     
175   if (!fPreCluster)
176   {
177     // we are done
178     return 0x0;
179   }
180     
181   WorkOnPreCluster();
182
183   // WorkOnPreCluster may have used only part of the pads, so we check that
184   // now, and let the unused pads be reused by the preclustering...
185   
186   Int_t mult = fPreCluster->Multiplicity();
187   for ( Int_t i = 0; i < mult; ++i )
188   {
189     AliMUONPad* pad = fPreCluster->Pad(i);
190     if ( !pad->IsUsed() )
191     {
192       fPreClusterFinder->UsePad(*pad);
193     }
194   }
195   
196   return NextCluster();
197 }
198
199 //_____________________________________________________________________________
200 Bool_t
201 AliMUONClusterFinderMLEM::WorkOnPreCluster()
202 {
203   /// Starting from a precluster, builds a pixel array, and then
204   /// extract clusters from this array
205   
206   //  AliCodeTimerAuto("",0)    
207
208   if (fDebug) {
209     cout << " *** Event # " << fEventNumber 
210          << " det. elem.: " << fDetElemId << endl;
211     for (Int_t j = 0; j < fPreCluster->Multiplicity(); ++j) {
212       AliMUONPad* pad = fPreCluster->Pad(j);
213       printf(" bbb %3d %1d %8.4f %8.4f %8.4f %8.4f %6.1f %3d %3d %2d %1d %1d \n",
214              j, pad->Cathode(), pad->Coord(0), pad->Coord(1), pad->DX()*2, pad->DY()*2,
215              pad->Charge(), pad->Ix(), pad->Iy(), pad->Status(), pad->IsReal(), pad->IsSaturated());
216     }
217   }
218
219   AliMUONCluster* cluster = CheckPrecluster(*fPreCluster);
220   if (!cluster) return kFALSE;
221
222   BuildPixArray(*cluster);
223   
224   if ( fPixArray->GetLast() < 0 )
225   {
226     AliDebug(1,"No pixel for the above cluster");
227     delete cluster;
228     return kFALSE;
229   }
230   
231   // Use MLEM for cluster finder
232   Int_t nMax = 1, localMax[100], maxPos[100];
233   Double_t maxVal[100];
234   
235   Int_t iSimple = 0, nInX = -1, nInY;
236   
237   PadsInXandY(*cluster,nInX, nInY);
238   
239   if (nInX < 4 && nInY < 4) 
240   {
241     iSimple = 1; // simple cluster
242   }
243   else 
244   {
245     nMax = FindLocalMaxima(fPixArray, localMax, maxVal); // for small clusters just to tag pixels
246     if (nMax > 1) {
247       if (cluster->Multiplicity() <= 50) nMax = 1; // for small clusters 
248       if (nMax > 1) TMath::Sort(nMax, maxVal, maxPos, kTRUE); // in descending order
249     }
250   }
251   
252   for (Int_t i = 0; i < nMax; ++i) 
253   {
254     if (nMax > 1) 
255     {
256       FindCluster(*cluster,localMax, maxPos[i]);
257     }
258
259     MainLoop(*cluster,iSimple);
260
261     if (i < nMax-1) 
262     {
263       Int_t mult = cluster->Multiplicity();
264       for (Int_t j = 0; j < mult; ++j) 
265       {
266         AliMUONPad* pad = cluster->Pad(j);
267         //if ( pad->Status() == 0 ) continue; // pad charge was not modified
268         if ( pad->Status() != fgkOver ) continue; // pad was not used
269         //pad->SetStatus(0);
270         pad->SetStatus(fgkZero);
271         pad->RevertCharge(); // use backup charge value
272       }
273     }
274   } // for (Int_t i=0; i<nMax;
275   delete fHistMlem;
276   delete fHistAnode;
277   fHistMlem = fHistAnode = 0x0;
278   delete cluster;
279   return kTRUE;
280 }
281
282 //_____________________________________________________________________________
283 Bool_t 
284 AliMUONClusterFinderMLEM::Overlap(const AliMUONPad& pad, const AliMUONPad& pixel)
285 {
286   /// Check if the pad and the pixel overlaps
287
288   // make a fake pad from the pixel
289   AliMUONPad tmp(pad.DetElemId(),pad.Cathode(),pad.Ix(),pad.Iy(),
290                  pixel.Coord(0),pixel.Coord(1),
291                  pixel.Size(0),pixel.Size(1),0);
292   
293   return AliMUONPad::AreOverlapping(pad,tmp,fgkDecreaseSize);
294 }
295
296 //_____________________________________________________________________________
297 AliMUONCluster* 
298 AliMUONClusterFinderMLEM::CheckPrecluster(const AliMUONCluster& origCluster)
299 {
300   /// Check precluster in order to attempt to simplify it (mostly for
301   /// two-cathode preclusters)
302     
303   AliCodeTimerAuto("",0)
304
305   // Disregard small clusters (leftovers from splitting or noise)
306   if ((origCluster.Multiplicity()==1 || origCluster.Multiplicity()==2) &&
307       origCluster.Charge(0)+origCluster.Charge(1) < fLowestClusterCharge )
308   { 
309     return 0x0;
310   }
311
312   AliMUONCluster* cluster = new AliMUONCluster(origCluster);
313
314   AliDebug(2,"Start of CheckPreCluster=");
315   //StdoutToAliDebug(2,cluster->Print("full"));
316
317   AliMUONCluster* rv(0x0);
318   
319   if (cluster->Multiplicity(0) && cluster->Multiplicity(1))
320   { 
321     rv = CheckPreclusterTwoCathodes(cluster);
322   }
323   else
324   {
325     rv = cluster;
326   }
327   return rv;
328 }
329
330 //_____________________________________________________________________________
331 AliMUONCluster*
332 AliMUONClusterFinderMLEM::CheckPreclusterTwoCathodes(AliMUONCluster* cluster)
333 {
334   /// Check two-cathode cluster
335   
336   Int_t npad = cluster->Multiplicity();
337   Int_t* flags = new Int_t[npad];
338   for (Int_t j = 0; j < npad; ++j) flags[j] = 0;
339   
340   // Check pad overlaps
341   for ( Int_t i = 0; i < npad; ++i) 
342   {
343     AliMUONPad* padi = cluster->Pad(i);
344     if ( padi->Cathode() != 0 ) continue;
345     for (Int_t j = i+1; j < npad; ++j) 
346     {
347       AliMUONPad* padj = cluster->Pad(j);
348       if ( padj->Cathode() != 1 ) continue;
349       if ( !AliMUONPad::AreOverlapping(*padi,*padj,fgkDecreaseSize) ) continue;
350       flags[i] = flags[j] = 1; // mark overlapped pads
351     } 
352   } 
353   
354   // Check if all pads overlap
355   Int_t nFlags=0;
356   for (Int_t i = 0; i < npad; ++i) 
357   {
358     if (!flags[i]) ++nFlags;
359   }
360   
361   if (nFlags > 0) 
362   {
363     // not all pads overlap.
364     if (fDebug) cout << " nFlags: " << nFlags << endl;
365     TObjArray toBeRemoved;
366     for (Int_t i = 0; i < npad; ++i) 
367     {
368       AliMUONPad* pad = cluster->Pad(i);
369       if (flags[i]) continue;
370       Int_t cath = pad->Cathode();
371       Int_t cath1 = TMath::Even(cath);
372       // Check for edge effect (missing pads on the _other_ cathode)
373       AliMpPad mpPad =
374       fkSegmentation[cath1]->PadByPosition(pad->Position().X(),
375                                            pad->Position().Y(),kFALSE);
376       if (!mpPad.IsValid()) continue;
377       if (nFlags == 1 && pad->Charge() < fLowestPadCharge) continue; 
378       AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
379                       fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),pad->Charge()));
380       toBeRemoved.AddLast(pad);
381       fPreCluster->Pad(i)->Release();
382     }
383     Int_t nRemove = toBeRemoved.GetEntriesFast();
384     for ( Int_t i = 0; i < nRemove; ++i )
385     {
386       cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
387     }
388   } 
389   
390   // Check correlations of cathode charges
391   if ( !cluster->IsSaturated() && cluster->ChargeAsymmetry() > 1 )
392   {
393     // big difference
394     Int_t cathode = cluster->MaxRawChargeCathode();
395     Int_t imin(-1);
396     Int_t imax(-1);
397     Double_t cmax(0);
398     Double_t cmin(1E9);
399     
400     // get min and max pad charges on the cathode opposite to the 
401     // max pad (given by MaxRawChargeCathode())
402     //
403     Int_t mult = cluster->Multiplicity();
404     for ( Int_t i = 0; i < mult; ++i )
405     {
406       AliMUONPad* pad = cluster->Pad(i);
407       if ( pad->Cathode() != cathode || !pad->IsReal() )
408       {
409         // only consider pads in the opposite cathode, and
410         // only consider real pads (i.e. exclude the virtual ones)
411         continue;
412       }
413       if ( pad->Charge() < cmin )
414       {
415         cmin = pad->Charge();
416         imin = i;
417         if (imax < 0) {
418           imax = imin;
419           cmax = cmin;
420         }
421       }
422       else if ( pad->Charge() > cmax )
423       {
424         cmax = pad->Charge();
425         imax = i;
426       }      
427     }
428     AliDebug(2,Form("Pad imin,imax %d,%d cmin,cmax %e,%e",
429                     imin,imax,cmin,cmax));
430     //
431     // arrange pads according to their distance to the max, normalized
432     // to the pad size
433     Double_t* dist = new Double_t[mult];
434     Double_t dxMin(1E9);
435     Double_t dyMin(1E9);
436     Double_t dmin(0);
437     
438     AliMUONPad* padmax = cluster->Pad(imax);
439     
440     for ( Int_t i = 0; i < mult; ++i )
441     {
442       dist[i] = 0.0;
443       if ( i == imax) continue;
444       AliMUONPad* pad = cluster->Pad(i);
445       if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
446       Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
447       Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
448       dist[i] = TMath::Sqrt(dx*dx+dy*dy);
449       if ( i == imin )
450       {
451         dmin = dist[i] + 1E-3; // distance to the pad with minimum charge
452         dxMin = dx;
453         dyMin = dy;
454       }      
455     }
456     
457     TMath::Sort(mult,dist,flags,kFALSE); // in ascending order
458     Double_t xmax(-1), distPrev(999);
459     TObjArray toBeRemoved;
460     
461     for ( Int_t i = 0; i < mult; ++i )
462     {
463       Int_t indx = flags[i];
464       AliMUONPad* pad = cluster->Pad(indx);
465       if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
466       if ( dist[indx] > dmin )
467       {
468         // farther than the minimum pad
469         Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
470         Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
471         dx *= dxMin;
472         dy *= dyMin;
473         if (dx >= 0 && dy >= 0) continue;
474         if (TMath::Abs(dx) > TMath::Abs(dy) && dx >= 0) continue;
475         if (TMath::Abs(dy) > TMath::Abs(dx) && dy >= 0) continue;        
476       }
477       if (dist[indx] > distPrev + 1) break; // overstepping empty pads
478       if ( pad->Charge() <= cmax || TMath::Abs(dist[indx]-xmax) < 1E-3 )
479       {
480         // release pad
481         if (TMath::Abs(dist[indx]-xmax) < 1.e-3) 
482         {
483           cmax = TMath::Max(pad->Charge(),cmax);
484         }
485         else
486         {
487           cmax = pad->Charge();
488         }
489         xmax = dist[indx];
490         distPrev = dist[indx];
491         AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
492                         fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),
493                         pad->Charge()));
494   
495         toBeRemoved.AddLast(pad);
496         fPreCluster->Pad(indx)->Release();
497       }
498     }
499     Int_t nRemove = toBeRemoved.GetEntriesFast();
500     for ( Int_t i = 0; i < nRemove; ++i )
501     {
502       cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
503     }    
504     delete[] dist;
505   }
506   
507   delete[] flags;
508   
509   AliDebug(2,"End of CheckPreClusterTwoCathodes=");
510   //StdoutToAliDebug(2,cluster->Print("full"));
511
512   return cluster;    
513 }
514
515 //_____________________________________________________________________________
516 void
517 AliMUONClusterFinderMLEM::CheckOverlaps()
518 {
519   /// For debug only : check if some pixels overlap...
520   
521   Int_t nPix = fPixArray->GetLast()+1;
522   Int_t dummy(0);
523   
524   for ( Int_t i = 0; i < nPix; ++i )
525   {
526     AliMUONPad* pixelI = Pixel(i);
527     AliMUONPad pi(dummy,dummy,dummy,dummy,
528                   pixelI->Coord(0),pixelI->Coord(1),
529                   pixelI->Size(0),pixelI->Size(1),0.0);
530     
531     for ( Int_t j = i+1; j < nPix; ++j )
532     {
533       AliMUONPad* pixelJ = Pixel(j);
534       AliMUONPad pj(dummy,dummy,dummy,dummy,
535                     pixelJ->Coord(0),pixelJ->Coord(1),
536                     pixelJ->Size(0),pixelJ->Size(1),0.0);  
537       AliMpArea area;
538       
539       if ( AliMUONPad::AreOverlapping(pi,pj,fgkDecreaseSize,area) )
540       {
541         AliInfo(Form("The following 2 pixels (%d and %d) overlap !",i,j));
542         /*
543         StdoutToAliInfo(pixelI->Print();
544                         cout << " Surface = " << pixelI->Size(0)*pixelI->Size(1)*4 << endl;
545                         pixelJ->Print();
546                         cout << " Surface = " << pixelJ->Size(0)*pixelJ->Size(1)*4 << endl;
547                         cout << " Area surface = " << area.GetDimensionX()*area.GetDimensionY()*4 << endl;
548                         cout << "-------" << endl;
549                         );
550         */        
551       }
552     }    
553   }
554 }
555
556 //_____________________________________________________________________________
557 void AliMUONClusterFinderMLEM::BuildPixArray(AliMUONCluster& cluster)
558 {
559   /// Build pixel array for MLEM method
560   
561   Int_t npad = cluster.Multiplicity();
562   if (npad<=0) 
563   {
564     AliWarning("Got no pad at all ?!");
565   }
566   
567   fPixArray->Delete();
568   BuildPixArrayOneCathode(cluster);
569   
570   Int_t nPix = fPixArray->GetLast()+1;
571   
572 //  AliDebug(2,Form("nPix after BuildPixArray=%d",nPix));
573   
574   if ( nPix > npad ) 
575   {
576 //    AliDebug(2,Form("Will trim number of pixels to number of pads"));
577     
578     // Too many pixels - sort and remove pixels with the lowest signal
579     fPixArray->Sort();
580     for ( Int_t i = npad; i < nPix; ++i ) 
581     {
582       RemovePixel(i);
583     }
584     fPixArray->Compress();
585   } // if (nPix > npad)
586
587 //  StdoutToAliDebug(2,cout << "End of BuildPixelArray:" << endl;
588 //                   fPixArray->Print(););
589   //CheckOverlaps();//FIXME : this is for debug only. Remove it.
590 }
591
592 //_____________________________________________________________________________
593 void AliMUONClusterFinderMLEM::BuildPixArrayOneCathode(AliMUONCluster& cluster)
594 {
595   /// Build the pixel array
596
597 //  AliDebug(2,Form("cluster.Multiplicity=%d",cluster.Multiplicity()));
598
599   TVector2 dim = cluster.MinPadDimensions (-1, kFALSE);
600   Double_t width[2] = {dim.X(), dim.Y()}, xy0[2]={99999,99999};
601   Int_t found[2] = {0,0}, mult = cluster.Multiplicity();
602
603   for ( Int_t i = 0; i < mult; ++i) {
604     AliMUONPad* pad = cluster.Pad(i);
605     for (Int_t j = 0; j < 2; ++j) {
606       if (found[j] == 0 && TMath::Abs(pad->Size(j)-width[j]) < fgkDistancePrecision) { 
607         xy0[j] = pad->Coord(j);
608         found[j] = 1;
609       }
610     }
611     if (found[0] && found[1]) break;
612   }
613
614   Double_t min[2], max[2];
615   Int_t cath0 = 0, cath1 = 1;
616   if (cluster.Multiplicity(0) == 0) cath0 = 1;
617   else if (cluster.Multiplicity(1) == 0) cath1 = 0;
618
619
620   Double_t leftDownX, leftDownY;
621   cluster.Area(cath0).LeftDownCorner(leftDownX, leftDownY);
622   Double_t rightUpX, rightUpY;
623   cluster.Area(cath0).RightUpCorner(rightUpX, rightUpY);
624   min[0] = leftDownX;
625   min[1] = leftDownY;
626   max[0] = rightUpX;
627   max[1] = rightUpY;;
628   if (cath1 != cath0) {
629     cluster.Area(cath1).LeftDownCorner(leftDownX, leftDownY);
630     cluster.Area(cath1).RightUpCorner(rightUpX, rightUpY);
631     min[0] = TMath::Max (min[0], leftDownX);
632     min[1] = TMath::Max (min[1], leftDownY);
633     max[0] = TMath::Min (max[0], rightUpX);
634     max[1] = TMath::Min (max[1], rightUpY);
635   }
636
637   // Adjust limits
638   //width[0] /= 2; width[1] /= 2; // just for check
639   Int_t nbins[2]={0,0};
640   for (Int_t i = 0; i < 2; ++i) {
641     Double_t dist = (min[i] - xy0[i]) / width[i] / 2;
642     if (TMath::Abs(dist) < 1.e-6) dist = -1.e-6;
643     min[i] = xy0[i] + (TMath::Nint(dist-TMath::Sign(1.e-6,dist)) 
644                        + TMath::Sign(0.5,dist)) * width[i] * 2;
645     nbins[i] = TMath::Nint ((max[i] - min[i]) / width[i] / 2);
646     if (nbins[i] == 0) ++nbins[i];
647     max[i] = min[i] + nbins[i] * width[i] * 2;
648     //cout << dist << " " << min[i] << " " << max[i] << " " << nbins[i] << endl;
649   }
650
651   // Book histogram
652   TH2D *hist1 = new TH2D ("Grid", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
653   TH2D *hist2 = new TH2D ("Entries", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
654   TAxis *xaxis = hist1->GetXaxis();
655   TAxis *yaxis = hist1->GetYaxis();
656
657   // Fill histogram
658   for ( Int_t i = 0; i < mult; ++i) {
659     AliMUONPad* pad = cluster.Pad(i);
660     Int_t ix0 = xaxis->FindBin(pad->X());
661     Int_t iy0 = yaxis->FindBin(pad->Y());
662     PadOverHist(0, ix0, iy0, pad, hist1, hist2);
663   }
664
665   // Store pixels
666   for (Int_t i = 1; i <= nbins[0]; ++i) {
667     Double_t x = xaxis->GetBinCenter(i);
668     for (Int_t j = 1; j <= nbins[1]; ++j) {
669       if (hist2->GetCellContent(i,j) < 0.1) continue;
670       //if (hist2->GetCellContent(i,j) < 1.1 && cluster.Multiplicity(0) && 
671       //  cluster.Multiplicity(1)) continue;
672       if (cath0 != cath1) {
673         // Two-sided cluster
674         Double_t cont = hist2->GetCellContent(i,j);
675         if (cont < 999.) continue;
676         if (cont-Int_t(cont/1000.)*1000. < 0.5) continue;
677       }
678       Double_t y = yaxis->GetBinCenter(j);
679       Double_t charge = hist1->GetCellContent(i,j);
680       AliMUONPad* pixPtr = new AliMUONPad(x, y, width[0], width[1], charge);
681       fPixArray->Add(pixPtr);
682     }  
683   }
684   //*
685   if (fPixArray->GetEntriesFast() == 1) {
686     // Split pixel into 2
687     AliMUONPad* pixPtr = static_cast<AliMUONPad*> (fPixArray->UncheckedAt(0));
688     pixPtr->SetSize(0,width[0]/2.);
689     pixPtr->Shift(0,-width[0]/4.);
690     pixPtr = new AliMUONPad(pixPtr->X()+width[0], pixPtr->Y(), width[0]/2., width[1], pixPtr->Charge());
691     fPixArray->Add(pixPtr);
692   }
693   //*/
694   //fPixArray->Print();
695   delete hist1;
696   delete hist2;
697 }
698
699 //_____________________________________________________________________________
700 void AliMUONClusterFinderMLEM::PadOverHist(Int_t idir, Int_t ix0, Int_t iy0, AliMUONPad *pad,
701                                            TH2D *hist1, TH2D *hist2)
702 {
703   /// "Span" pad over histogram in the direction idir
704
705   TAxis *axis = idir == 0 ? hist1->GetXaxis() : hist1->GetYaxis();
706   Int_t nbins = axis->GetNbins(), cath = pad->Cathode();
707   Double_t bin = axis->GetBinWidth(1), amask = TMath::Power(1000.,cath*1.);
708
709   Int_t nbinPad = (Int_t)(pad->Size(idir)/bin*2+fgkDistancePrecision) + 1; // number of bins covered by pad
710
711   for (Int_t i = 0; i < nbinPad; ++i) {
712     Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
713     if (ixy > nbins) break;
714     Double_t lowEdge = axis->GetBinLowEdge(ixy);
715     if (lowEdge + fgkDistancePrecision > pad->Coord(idir) + pad->Size(idir)) break;
716     if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
717     else {
718       // Fill histogram
719       Double_t cont = pad->Charge();
720       if (hist2->GetCellContent(ix0, ixy) > 0.1) 
721         cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont);
722       hist1->SetCellContent(ix0, ixy, cont);
723       //hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+1);
724       hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
725     }
726   }
727
728   for (Int_t i = -1; i > -nbinPad; --i) {
729     Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
730     if (ixy < 1) break;
731     Double_t upEdge = axis->GetBinUpEdge(ixy);
732     if (upEdge - fgkDistancePrecision < pad->Coord(idir) - pad->Size(idir)) break;
733     if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
734     else {
735       // Fill histogram
736       Double_t cont = pad->Charge();
737       if (hist2->GetCellContent(ix0, ixy) > 0.1) 
738         cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont);
739       hist1->SetCellContent(ix0, ixy, cont);
740       //hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+1);
741       hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
742     }
743   }
744 }
745
746 //_____________________________________________________________________________
747 void
748 AliMUONClusterFinderMLEM::Plot(const char* basename)
749 {
750   /// Make a plot and save it as png
751   
752   return; //AZ
753   if (!fPlot) return;
754   
755   TCanvas* c = new TCanvas("MLEM","MLEM",800,600);
756   c->Draw();
757   Draw();
758   c->Modified();
759   c->Update();
760   c->Print(Form("%s.EVT%d.DE%d.CLU%d.png",basename,fEventNumber,
761                 fDetElemId,fClusterNumber));
762 }
763
764 //_____________________________________________________________________________
765 void
766 AliMUONClusterFinderMLEM::ComputeCoefficients(AliMUONCluster& cluster,
767                                               Double_t* coef,
768                                               Double_t* probi)
769 {
770   /// Compute coefficients needed for MLEM algorithm
771   
772   Int_t npadTot = cluster.Multiplicity();
773   Int_t nPix = fPixArray->GetLast()+1;
774   
775   //memset(probi,0,nPix*sizeof(Double_t));
776   for (Int_t j = 0; j < npadTot*nPix; ++j) coef[j] = 0.;
777   for (Int_t j = 0; j < nPix; ++j) probi[j] = 0.;
778
779   Int_t mult = cluster.Multiplicity();
780   for ( Int_t j = 0; j < mult; ++j ) 
781   {
782     AliMUONPad* pad = cluster.Pad(j);
783     Int_t indx = j*nPix;
784   
785     for ( Int_t ipix = 0; ipix < nPix; ++ipix ) 
786     {
787       Int_t indx1 = indx + ipix;
788       //if (pad->Status() < 0) 
789       if (pad->Status() != fgkZero) 
790       {   
791         coef[indx1] = 0; 
792         continue; 
793       }
794       AliMUONPad* pixPtr = Pixel(ipix);
795       // coef is the charge (given by Mathieson integral) on pad, assuming
796       // the Mathieson is center at pixel.
797       coef[indx1] = fSplitter->ChargeIntegration(pixPtr->Coord(0), pixPtr->Coord(1), *pad);  
798 //      AliDebug(2,Form("pad=(%d,%d,%e,%e,%e,%e) pix=(%e,%e,%e,%e) coef %e",
799 //                      pad->Ix(),pad->Iy(),
800 //                      pad->X(),pad->Y(),
801 //                      pad->DX(),pad->DY(),
802 //                      pixPtr->Coord(0),pixPtr->Coord(1), 
803 //                      pixPtr->Size(0),pixPtr->Size(1),
804 //                      coef[indx1]));
805       
806       probi[ipix] += coef[indx1];
807     } 
808   } 
809 }
810
811 //_____________________________________________________________________________
812 Bool_t AliMUONClusterFinderMLEM::MainLoop(AliMUONCluster& cluster, Int_t iSimple)
813 {
814   /// Repeat MLEM algorithm until pixel size becomes sufficiently small
815   
816   //  AliCodeTimerAuto("",0)
817   
818   Int_t nPix = fPixArray->GetLast()+1;
819
820   AliDebug(2,Form("nPix=%d iSimple=%d, precluster=",nPix,iSimple));
821   //StdoutToAliDebug(2,cluster.Print("full"););
822
823   if ( nPix < 0 )
824   {
825     AliDebug(1,"No pixels, why am I here then ?");
826   }
827   
828   AddVirtualPad(cluster); // add virtual pads if necessary
829   
830   Int_t npadTot = cluster.Multiplicity();
831   Int_t npadOK = 0;
832   for (Int_t i = 0; i < npadTot; ++i) 
833   {
834     //if (cluster.Pad(i)->Status() == 0) ++npadOK;
835     if (cluster.Pad(i)->Status() == fgkZero) ++npadOK;
836   }
837
838   Double_t* coef(0x0);
839   Double_t* probi(0x0);
840   Int_t lc(0); // loop counter
841   
842   //Plot("mlem.start");
843   AliMUONPad* pixPtr = Pixel(0);
844   Double_t xylim[4] = {pixPtr->X(), -pixPtr->X(), pixPtr->Y(), -pixPtr->Y()};
845
846   while (1) 
847   {
848     ++lc;
849     delete fHistMlem;
850     fHistMlem = 0x0;
851     
852     AliDebug(2,Form("lc %d nPix %d(%d) npadTot %d npadOK %d",lc,nPix,fPixArray->GetLast()+1,npadTot,npadOK));
853     AliDebug(2,Form("EVT%d PixArray=",fEventNumber));
854     //StdoutToAliDebug(2,fPixArray->Print("full"));
855         
856     coef = new Double_t [npadTot*nPix];
857     probi = new Double_t [nPix];
858
859     // Calculate coefficients and pixel visibilities
860     ComputeCoefficients(cluster,coef,probi);
861
862     for (Int_t ipix = 0; ipix < nPix; ++ipix) 
863     {
864       if (probi[ipix] < 0.01) 
865       {
866         AliMUONPad* pixel = Pixel(ipix);
867         AliDebug(2,Form("Setting the following pixel to invisible as its probi<0.01:"));
868         //StdoutToAliDebug(2,cout << Form(" -- ipix %3d --- "); pixel->Print(););
869         pixel->SetCharge(0); // "invisible" pixel
870       }
871     }
872     
873     // MLEM algorithm
874     Mlem(cluster,coef, probi, 15);
875
876     // Find histogram limits for the 1'st pass only - for others computed below
877     if (lc == 1) {
878       for ( Int_t ipix = 1; ipix < nPix; ++ipix ) 
879         {
880           pixPtr = Pixel(ipix);
881           for ( Int_t i = 0; i < 2; ++i ) 
882             {
883               Int_t indx = i * 2;
884               if (pixPtr->Coord(i) < xylim[indx]) xylim[indx] = pixPtr->Coord(i); 
885               else if (-pixPtr->Coord(i) < xylim[indx+1]) xylim[indx+1] = -pixPtr->Coord(i); 
886             }
887         }
888     } else pixPtr = Pixel(0);
889
890     for (Int_t i = 0; i < 4; i++) 
891     {
892       xylim[i] -= pixPtr->Size(i/2); 
893     }
894     
895     Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
896     Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
897
898     //StdoutToAliDebug(2,cout << "pixel used for nx,ny computation : "; pixPtr->Print(););
899     AliDebug(2,Form("lc %d pixPtr size = %e,%e nx,ny=%d,%d xylim=%e,%e,%e,%e",
900                     lc,pixPtr->Size(0),pixPtr->Size(1),nx,ny,
901                     xylim[0],-xylim[1],xylim[2],-xylim[3]
902                     ));
903     
904     AliDebug(2,Form("LowestPadCharge=%e",fLowestPadCharge));
905     
906     delete fHistMlem;
907     
908     fHistMlem = new TH2D("mlem","mlem",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
909
910     for (Int_t ipix = 0; ipix < nPix; ++ipix) 
911     {
912       AliMUONPad* pixPtr2 = Pixel(ipix);
913       fHistMlem->Fill(pixPtr2->Coord(0),pixPtr2->Coord(1),pixPtr2->Charge());
914     }
915
916     // Check if the total charge of pixels is too low
917     Double_t qTot = 0;
918     for ( Int_t i = 0; i < nPix; ++i) 
919     {
920       qTot += Pixel(i)->Charge();
921     }
922     
923     if ( qTot < 1.e-4 || ( npadOK < 3 && qTot < fLowestClusterCharge ) )
924     {
925       AliDebug(1,Form("Deleting the above cluster (charge %e too low, npadOK=%d)",qTot,npadOK));
926       delete [] coef; 
927       delete [] probi; 
928       fPixArray->Delete(); 
929       for ( Int_t i = 0; i < npadTot; ++i) 
930       {
931         AliMUONPad* pad = cluster.Pad(i);
932         //if ( pad->Status() == 0) pad->SetStatus(-1);
933         if ( pad->Status() == fgkZero) pad->SetStatus(fgkOver);
934       }
935       return kFALSE; 
936     }
937
938     if (iSimple) 
939     {
940       // Simple cluster - skip further passes thru EM-procedure
941       Simple(cluster);
942       delete [] coef; 
943       delete [] probi; 
944       fPixArray->Delete(); 
945       return kTRUE;
946     }
947
948     // Calculate position of the center-of-gravity around the maximum pixel
949     Double_t xyCOG[2];
950     FindCOG(xyCOG);
951
952     if (TMath::Min(pixPtr->Size(0),pixPtr->Size(1)) < 0.07 && 
953         pixPtr->Size(0) > pixPtr->Size(1)) break;
954
955     // Sort pixels according to the charge
956     MaskPeaks(1); // mask local maxima
957     fPixArray->Sort();
958     MaskPeaks(0); // unmask local maxima
959     Double_t pixMin = 0.01*Pixel(0)->Charge();
960     pixMin = TMath::Min(pixMin,100*fLowestPixelCharge);
961
962     // Decrease pixel size and shift pixels to make them centered at 
963     // the maximum one
964     Int_t indx = (pixPtr->Size(0)>pixPtr->Size(1)) ? 0 : 1;
965     Int_t ix(1);
966     Double_t width = 0;
967     Double_t shift[2] = { 0.0, 0.0 };
968     for (Int_t i = 0; i < 4; ++i) xylim[i] = 999;
969     Int_t nPix1 = nPix; 
970     nPix = 0;
971     for (Int_t ipix = 0; ipix < nPix1; ++ipix) 
972     {
973       AliMUONPad* pixPtr2 = Pixel(ipix);
974       if ( nPix >= npadOK  // too many pixels already
975            ||
976            ((pixPtr2->Charge() < pixMin) && (pixPtr2->Status() != fgkMustKeep)) // too low charge
977            ) 
978       { 
979         RemovePixel(ipix);
980         continue;
981       }
982       for (Int_t i = 0; i < 2; ++i) 
983       {
984         if (!i) 
985         {
986           pixPtr2->SetCharge(fLowestPadCharge);
987           pixPtr2->SetSize(indx, pixPtr2->Size(indx)/2);
988           width = -pixPtr2->Size(indx);
989           pixPtr2->Shift(indx, width);
990           // Shift pixel position
991           if (ix) 
992           {
993             ix = 0;
994             for (Int_t j = 0; j < 2; ++j) 
995             {
996               shift[j] = pixPtr2->Coord(j) - xyCOG[j];
997               shift[j] -= ((Int_t)(shift[j]/pixPtr2->Size(j)/2))*pixPtr2->Size(j)*2;
998             }
999           } // if (ix)
1000           pixPtr2->Shift(0, -shift[0]);
1001           pixPtr2->Shift(1, -shift[1]);
1002           ++nPix;
1003         } 
1004         else if (nPix < npadOK)
1005         {
1006           pixPtr2 = new AliMUONPad(*pixPtr2);
1007           pixPtr2->Shift(indx, -2*width);
1008           pixPtr2->SetStatus(fgkZero);
1009           fPixArray->Add(pixPtr2);
1010           ++nPix;
1011         } 
1012         else continue; // skip adjustment of histo limits
1013         for (Int_t j = 0; j < 4; ++j) 
1014         {
1015           xylim[j] = TMath::Min (xylim[j], (j%2 ? -1 : 1)*pixPtr2->Coord(j/2));
1016         }
1017       } // for (Int_t i=0; i<2;
1018     } // for (Int_t ipix=0;
1019     
1020     fPixArray->Compress();
1021
1022     AliDebug(2,Form("After shift:"));
1023     //StdoutToAliDebug(2,fPixArray->Print("","full"););
1024     //Plot(Form("mlem.lc%d",lc+1));
1025     
1026     AliDebug(2,Form(" xyCOG=%9.6f %9.6f xylim=%9.6f,%9.6f,%9.6f,%9.6f",
1027                     xyCOG[0],xyCOG[1],
1028                     xylim[0],xylim[1],
1029                     xylim[2],xylim[3]));
1030
1031     if (nPix < npadOK)
1032     {
1033       AliMUONPad* pixPtr2 = Pixel(0);
1034       // add pixels if the maximum is at the limit of pixel area:
1035       // start from Y-direction
1036       Int_t j = 0;
1037       for (Int_t i = 3; i > -1; --i) 
1038       {
1039         if (nPix < npadOK && 
1040             TMath::Abs((i%2 ? -1 : 1)*xylim[i]-xyCOG[i/2]) < pixPtr2->Size(i/2)) 
1041         {
1042           //AliMUONPad* p = static_cast<AliMUONPad*>(pixPtr->Clone());
1043           AliMUONPad* p = new AliMUONPad(*pixPtr2);
1044           p->SetCoord(i/2, xyCOG[i/2]+(i%2 ? 2:-2)*pixPtr2->Size(i/2));
1045           xylim[i] = p->Coord(i/2) * (i%2 ? -1 : 1); // update histo limits
1046           j = TMath::Even (i/2);
1047           p->SetCoord(j, xyCOG[j]);
1048           AliDebug(2,Form("Adding pixel on the edge (i=%d) ",i));
1049           //StdoutToAliDebug(2,cout << " ---- "; 
1050           //               p->Print("corners"););
1051           fPixArray->Add(p);
1052           ++nPix;
1053         }
1054       }
1055     } 
1056     nPix = fPixArray->GetEntriesFast();
1057     delete [] coef; 
1058     delete [] probi; 
1059   } // while (1)
1060
1061   AliDebug(2,Form("At the end of while loop nPix=%d : ",fPixArray->GetLast()+1));
1062   //StdoutToAliDebug(2,fPixArray->Print("","full"););
1063
1064   // remove pixels with low signal or low visibility
1065   // Cuts are empirical !!!
1066   Double_t thresh = TMath::Max (fHistMlem->GetMaximum()/100.,2.0*fLowestPixelCharge);
1067   thresh = TMath::Min (thresh,100.0*fLowestPixelCharge);
1068   Double_t charge = 0;
1069
1070   // Mark pixels which should be removed
1071   for (Int_t i = 0; i < nPix; ++i) 
1072   {
1073     AliMUONPad* pixPtr2 = Pixel(i);
1074     charge = pixPtr2->Charge();
1075     if (charge < thresh) 
1076     {
1077       pixPtr2->SetCharge(-charge);
1078     }
1079   }
1080
1081   // Move charge of removed pixels to their nearest neighbour (to keep total charge the same)
1082   Int_t near = 0;
1083   for (Int_t i = 0; i < nPix; ++i) 
1084   {
1085     AliMUONPad* pixPtr2 = Pixel(i);
1086     charge = pixPtr2->Charge();
1087     if (charge > 0) continue;
1088     near = FindNearest(pixPtr2);
1089     pixPtr2->SetCharge(0);
1090     probi[i] = 0; // make it "invisible"
1091     AliMUONPad* pnear = Pixel(near);
1092     pnear->SetCharge(pnear->Charge() + (-charge));
1093   }
1094   Mlem(cluster,coef,probi,2);
1095   
1096   AliDebug(2,Form("Before splitting nPix=%d EVT %d DE %d",fPixArray->GetLast()+1,fEventNumber,fDetElemId));
1097   //StdoutToAliDebug(2,fPixArray->Print("","full"););
1098   //Plot("mlem.beforesplit");
1099   
1100   // Update histogram
1101   for (Int_t i = 0; i < nPix; ++i) 
1102   {
1103     AliMUONPad* pixPtr2 = Pixel(i);
1104     Int_t ix = fHistMlem->GetXaxis()->FindBin(pixPtr2->Coord(0));
1105     Int_t iy = fHistMlem->GetYaxis()->FindBin(pixPtr2->Coord(1));
1106     fHistMlem->SetBinContent(ix, iy, pixPtr2->Charge());
1107   }
1108
1109   // Try to split into clusters
1110   Bool_t ok = kTRUE;
1111   if (fHistMlem->GetSum() < 2.0*fLowestPixelCharge) 
1112   {
1113     ok = kFALSE;
1114   }
1115   else 
1116   {
1117     fSplitter->Split(cluster,fHistMlem,coef,fClusterList);
1118   }
1119   
1120   delete [] coef; 
1121   delete [] probi; 
1122   fPixArray->Delete(); 
1123   
1124   return ok;
1125 }
1126
1127 //_____________________________________________________________________________
1128 void AliMUONClusterFinderMLEM::MaskPeaks(Int_t mask)
1129 {
1130   /// Mask/unmask pixels corresponding to local maxima (add/subtract 10000 to their charge
1131   /// - to avoid loosing low charge pixels after sorting)
1132
1133   for (Int_t i = 0; i < fPixArray->GetEntriesFast(); ++i) {
1134     AliMUONPad* pix = Pixel(i);
1135     if (pix->Status() == fgkMustKeep) {
1136       if (mask == 1) pix->SetCharge(pix->Charge()+10000.);
1137       else pix->SetCharge(pix->Charge()-10000.);
1138     }
1139   }
1140 }
1141
1142 //_____________________________________________________________________________
1143 void AliMUONClusterFinderMLEM::Mlem(AliMUONCluster& cluster, 
1144                                     const Double_t* coef, Double_t* probi, 
1145                                     Int_t nIter)
1146 {
1147   /// Use MLEM to find pixel charges
1148   
1149   Int_t nPix = fPixArray->GetEntriesFast();
1150
1151   Int_t npad = cluster.Multiplicity();
1152
1153   Double_t* probi1 = new Double_t[nPix];
1154   Double_t probMax = TMath::MaxElement(nPix,probi);
1155   
1156   for (Int_t iter = 0; iter < nIter; ++iter) 
1157   {
1158     // Do iterations
1159     for (Int_t ipix = 0; ipix < nPix; ++ipix) 
1160     {
1161       Pixel(ipix)->SetChargeBackup(0);
1162       // Correct each pixel
1163       probi1[ipix] = 0;
1164       if (probi[ipix] < 0.01) continue; // skip "invisible" pixel
1165       Double_t sum = 0;
1166       probi1[ipix] = probMax;
1167       for (Int_t j = 0; j < npad; ++j) 
1168       {
1169         AliMUONPad* pad = cluster.Pad(j);
1170         //if (pad->Status() < 0) continue; 
1171         if (pad->Status() != fgkZero) continue; 
1172         Double_t sum1 = 0;
1173         Int_t indx1 = j*nPix;
1174         Int_t indx = indx1 + ipix;
1175         // Calculate expectation
1176         for (Int_t i = 0; i < nPix; ++i) 
1177         {
1178           sum1 += Pixel(i)->Charge()*coef[indx1+i];
1179           //cout << i << " " << Pixel(i)->Charge() << " " << coef[indx1+i] << endl;
1180         } 
1181         if ( pad->IsSaturated() && sum1 > pad->Charge() ) 
1182         { 
1183           // correct for pad charge overflows
1184           probi1[ipix] -= coef[indx]; 
1185           continue; 
1186           //sum1 = pad->Charge();
1187         } 
1188
1189         if (sum1 > 1.e-6) sum += pad->Charge()*coef[indx]/sum1;
1190       } // for (Int_t j=0;
1191       AliMUONPad* pixPtr = Pixel(ipix);
1192       if (probi1[ipix] > 1.e-6) 
1193       {
1194         //AZ pixPtr->SetCharge(pixPtr->Charge()*sum/probi1[ipix]);
1195         pixPtr->SetChargeBackup(pixPtr->Charge()*sum/probi1[ipix]);
1196       }
1197       //cout << " xxx " << ipix << " " << pixPtr->Charge() << " " << pixPtr->ChargeBackup() << " " << sum << " " << probi1[ipix] << endl;
1198     } // for (Int_t ipix=0;
1199     Double_t qTot = 0;
1200     for (Int_t i = 0; i < nPix; ++i) {
1201       AliMUONPad* pixPtr = Pixel(i);
1202       pixPtr->RevertCharge();
1203       qTot += pixPtr->Charge();
1204     }
1205     if (qTot < 1.e-6) {
1206       // Can happen in clusters with large number of overflows - speeding up 
1207       delete [] probi1;
1208       return;
1209     }
1210   } // for (Int_t iter=0;
1211   delete [] probi1;
1212 }
1213
1214 //_____________________________________________________________________________
1215 void AliMUONClusterFinderMLEM::FindCOG(Double_t *xyc)
1216 {
1217   /// Calculate position of the center-of-gravity around the maximum pixel
1218
1219   Int_t ixmax, iymax, ix, nsumx=0, nsumy=0, nsum=0;
1220   Int_t i1 = -9, j1 = -9;
1221   fHistMlem->GetMaximumBin(ixmax,iymax,ix);
1222   Int_t nx = fHistMlem->GetNbinsX();
1223   Int_t ny = fHistMlem->GetNbinsY();
1224   Double_t thresh = fHistMlem->GetMaximum()/10;
1225   Double_t x, y, cont, xq=0, yq=0, qq=0;
1226   
1227   Int_t ie = TMath::Min(ny,iymax+1), je = TMath::Min(nx,ixmax+1);
1228   for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
1229     y = fHistMlem->GetYaxis()->GetBinCenter(i);
1230     for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
1231       cont = fHistMlem->GetCellContent(j,i);
1232       if (cont < thresh) continue;
1233       if (i != i1) {i1 = i; nsumy++;}
1234       if (j != j1) {j1 = j; nsumx++;}
1235       x = fHistMlem->GetXaxis()->GetBinCenter(j);
1236       xq += x*cont;
1237       yq += y*cont;
1238       qq += cont;
1239       nsum++;
1240     }
1241   }
1242   
1243   Double_t cmax = 0;
1244   Int_t i2 = 0, j2 = 0;
1245   x = y = 0;
1246   if (nsumy == 1) {
1247     // one bin in Y - add one more (with the largest signal)
1248     for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
1249       if (i == iymax) continue;
1250       for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
1251         cont = fHistMlem->GetCellContent(j,i);
1252         if (cont > cmax) {
1253           cmax = cont;
1254           x = fHistMlem->GetXaxis()->GetBinCenter(j);
1255           y = fHistMlem->GetYaxis()->GetBinCenter(i);
1256           i2 = i;
1257           j2 = j;
1258         }
1259       }
1260     }
1261     xq += x*cmax;
1262     yq += y*cmax;
1263     qq += cmax;
1264     if (i2 != i1) nsumy++;
1265     if (j2 != j1) nsumx++;
1266     nsum++;
1267   } // if (nsumy == 1)
1268   
1269   if (nsumx == 1) {
1270     // one bin in X - add one more (with the largest signal)
1271     cmax = x = y = 0;
1272     for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
1273       if (j == ixmax) continue;
1274       for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
1275         cont = fHistMlem->GetCellContent(j,i);
1276         if (cont > cmax) {
1277           cmax = cont;
1278           x = fHistMlem->GetXaxis()->GetBinCenter(j);
1279           y = fHistMlem->GetYaxis()->GetBinCenter(i);
1280           i2 = i;
1281           j2 = j;
1282         }
1283       }
1284     }
1285     xq += x*cmax;
1286     yq += y*cmax;
1287     qq += cmax;
1288     if (i2 != i1) nsumy++;
1289     if (j2 != j1) nsumx++;
1290     nsum++;
1291   } // if (nsumx == 1)
1292   
1293   xyc[0] = xq/qq; xyc[1] = yq/qq;
1294   AliDebug(2,Form("x,y COG = %e,%e",xyc[0],xyc[1]));
1295 }
1296
1297 //_____________________________________________________________________________
1298 Int_t AliMUONClusterFinderMLEM::FindNearest(const AliMUONPad *pixPtr0)
1299 {
1300 /// Find the pixel nearest to the given one
1301 /// (algorithm may be not very efficient)
1302
1303   Int_t nPix = fPixArray->GetEntriesFast(), imin = 0;
1304   Double_t rmin = 99999, dx = 0, dy = 0, r = 0;
1305   Double_t xc = pixPtr0->Coord(0), yc = pixPtr0->Coord(1);
1306   AliMUONPad *pixPtr;
1307
1308   for (Int_t i = 0; i < nPix; ++i) {
1309     pixPtr = (AliMUONPad*) fPixArray->UncheckedAt(i);
1310     if (pixPtr == pixPtr0 || pixPtr->Charge() < fLowestPixelCharge) continue;
1311     dx = (xc - pixPtr->Coord(0)) / pixPtr->Size(0);
1312     dy = (yc - pixPtr->Coord(1)) / pixPtr->Size(1);
1313     r = dx *dx + dy * dy;
1314     if (r < rmin) { rmin = r; imin = i; }
1315   }
1316   return imin;
1317 }
1318
1319 //_____________________________________________________________________________
1320 void
1321 AliMUONClusterFinderMLEM::Paint(Option_t*)
1322 {
1323   /// Paint cluster and pixels
1324   
1325   AliMpArea area(fPreCluster->Area());
1326   
1327   gPad->Range(area.LeftBorder(),area.DownBorder(),area.RightBorder(),area.UpBorder());
1328
1329   gVirtualX->SetFillStyle(1001);
1330   gVirtualX->SetFillColor(3);    
1331   gVirtualX->SetLineColor(3);
1332   
1333   Double_t s(1.0);
1334   
1335   for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1336   {
1337     AliMUONPad* pixel = Pixel(i);
1338
1339     gPad->PaintBox(pixel->Coord(0)-pixel->Size(0)*s,
1340                    pixel->Coord(1)-pixel->Size(1)*s,
1341                    pixel->Coord(0)+pixel->Size(0)*s,
1342                    pixel->Coord(1)+pixel->Size(1)*s);
1343
1344 //    for ( Int_t sign = -1; sign < 2; sign +=2 )
1345 //    {
1346 //      gPad->PaintLine(pixel->Coord(0) - pixel->Size(0),
1347 //                      pixel->Coord(1) + sign*pixel->Size(1),
1348 //                      pixel->Coord(0) + pixel->Size(0),
1349 //                      pixel->Coord(1) - sign*pixel->Size(1)
1350 //                    );
1351 //    }
1352   }      
1353
1354
1355   gVirtualX->SetFillStyle(0);
1356   
1357   fPreCluster->Paint();
1358
1359   gVirtualX->SetLineColor(1);
1360   gVirtualX->SetLineWidth(2);
1361   gVirtualX->SetFillStyle(0);
1362   gVirtualX->SetTextColor(1);
1363   gVirtualX->SetTextAlign(22);
1364   
1365   for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1366   {
1367     AliMUONPad* pixel = Pixel(i);
1368     gPad->PaintBox(pixel->Coord(0)-pixel->Size(0),
1369                    pixel->Coord(1)-pixel->Size(1),
1370                    pixel->Coord(0)+pixel->Size(0),
1371                    pixel->Coord(1)+pixel->Size(1));    
1372     gVirtualX->SetTextSize(pixel->Size(0)*60);
1373
1374     gPad->PaintText(pixel->Coord(0),pixel->Coord(1),Form("%d",(Int_t)(pixel->Charge())));
1375   }  
1376 }
1377
1378 //_____________________________________________________________________________
1379 Int_t AliMUONClusterFinderMLEM::FindLocalMaxima(TObjArray *pixArray, Int_t *localMax, Double_t *maxVal)
1380 {
1381 /// Find local maxima in pixel space for large preclusters in order to
1382 /// try to split them into smaller pieces (to speed up the MLEM procedure)
1383 /// or to find additional fitting seeds if clusters were not completely resolved  
1384
1385   AliDebug(1,Form("nPix=%d",pixArray->GetLast()+1));
1386
1387   Double_t xylim[4] = {999, 999, 999, 999};
1388
1389   Int_t nPix = pixArray->GetEntriesFast();
1390   AliMUONPad *pixPtr = 0;
1391   for (Int_t ipix = 0; ipix < nPix; ++ipix) {
1392     pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
1393     for (Int_t i = 0; i < 4; ++i) 
1394          xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
1395   }
1396   for (Int_t i = 0; i < 4; ++i) xylim[i] -= pixPtr->Size(i/2); 
1397
1398   Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
1399   Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
1400   if (pixArray == fPixArray) fHistAnode = new TH2D("anode","anode",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1401   else fHistAnode = new TH2D("anode1","anode1",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1402   for (Int_t ipix = 0; ipix < nPix; ++ipix) {
1403     pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
1404     fHistAnode->Fill(pixPtr->Coord(0), pixPtr->Coord(1), pixPtr->Charge());
1405   }
1406 //  if (fDraw && pixArray == fPixArray) fDraw->DrawHist("c2", hist);
1407
1408   Int_t nMax = 0, indx, nxy = ny * nx;
1409   Int_t *isLocalMax = new Int_t[nxy];
1410   for (Int_t i = 0; i < nxy; ++i) isLocalMax[i] = 0; 
1411
1412   for (Int_t i = 1; i <= ny; ++i) {
1413     indx = (i-1) * nx;
1414     for (Int_t j = 1; j <= nx; ++j) {
1415       if (fHistAnode->GetCellContent(j,i) < fLowestPixelCharge) continue;
1416       //if (isLocalMax[indx+j-1] < 0) continue;
1417       if (isLocalMax[indx+j-1] != 0) continue;
1418       FlagLocalMax(fHistAnode, i, j, isLocalMax);
1419     }
1420   }
1421
1422   for (Int_t i = 1; i <= ny; ++i) {
1423     indx = (i-1) * nx;
1424     for (Int_t j = 1; j <= nx; ++j) {
1425       if (isLocalMax[indx+j-1] > 0) { 
1426         localMax[nMax] = indx + j - 1; 
1427         maxVal[nMax++] = fHistAnode->GetCellContent(j,i);
1428         ((AliMUONPad*)fSplitter->BinToPix(fHistAnode, j, i))->SetStatus(fgkMustKeep);
1429         if (nMax > 99) break;
1430       }
1431     }
1432     if (nMax > 99) {
1433       AliError(" Too many local maxima !!!");
1434       break;
1435     }
1436   }
1437   if (fDebug) cout << " Local max: " << nMax << endl;
1438   delete [] isLocalMax; 
1439   return nMax;
1440 }
1441
1442 //_____________________________________________________________________________
1443 void AliMUONClusterFinderMLEM::FlagLocalMax(TH2D *hist, Int_t i, Int_t j, Int_t *isLocalMax)
1444 {
1445 /// Flag pixels (whether or not local maxima)
1446
1447   Int_t nx = hist->GetNbinsX();
1448   Int_t ny = hist->GetNbinsY();
1449   Int_t cont = TMath::Nint (hist->GetCellContent(j,i));
1450   Int_t cont1 = 0, indx = (i-1)*nx+j-1, indx1 = 0, indx2 = 0;
1451
1452   Int_t ie = i + 2, je = j + 2;
1453   for (Int_t i1 = i-1; i1 < ie; ++i1) {
1454     if (i1 < 1 || i1 > ny) continue;
1455     indx1 = (i1 - 1) * nx;
1456     for (Int_t j1 = j-1; j1 < je; ++j1) {
1457       if (j1 < 1 || j1 > nx) continue;
1458       if (i == i1 && j == j1) continue;
1459       indx2 = indx1 + j1 - 1;
1460       cont1 = TMath::Nint (hist->GetCellContent(j1,i1));
1461       if (cont < cont1) { isLocalMax[indx] = -1; return; }
1462       else if (cont > cont1) isLocalMax[indx2] = -1;
1463       else { // the same charge
1464         isLocalMax[indx] = 1; 
1465         if (isLocalMax[indx2] == 0) {
1466           FlagLocalMax(hist, i1, j1, isLocalMax);
1467           if (isLocalMax[indx2] < 0) { isLocalMax[indx] = -1; return; }
1468           else isLocalMax[indx2] = -1;
1469         }
1470       } 
1471     }
1472   }
1473   isLocalMax[indx] = 1; // local maximum
1474 }
1475
1476 //_____________________________________________________________________________
1477 void AliMUONClusterFinderMLEM::FindCluster(AliMUONCluster& cluster, 
1478                                            const Int_t *localMax, Int_t iMax)
1479 {
1480 /// Find pixel cluster around local maximum \a iMax and pick up pads
1481 /// overlapping with it
1482
1483   /* Just for check
1484   TCanvas* c = new TCanvas("Anode","Anode",800,600);
1485   c->cd();
1486   hist->Draw("lego1Fb"); // debug
1487   c->Update();
1488   Int_t tmp;
1489   cin >> tmp;
1490   */
1491   Int_t nx = fHistAnode->GetNbinsX();
1492   Int_t ny = fHistAnode->GetNbinsY();
1493   Int_t ic = localMax[iMax] / nx + 1;
1494   Int_t jc = localMax[iMax] % nx + 1;
1495   Int_t nxy = ny * nx;
1496   Bool_t *used = new Bool_t[nxy];
1497   for (Int_t i = 0; i < nxy; ++i) used[i] = kFALSE; 
1498
1499   // Drop all pixels from the array - pick up only the ones from the cluster
1500   fPixArray->Delete();
1501
1502   Double_t wx = fHistAnode->GetXaxis()->GetBinWidth(1)/2; 
1503   Double_t wy = fHistAnode->GetYaxis()->GetBinWidth(1)/2;  
1504   Double_t yc = fHistAnode->GetYaxis()->GetBinCenter(ic);
1505   Double_t xc = fHistAnode->GetXaxis()->GetBinCenter(jc);
1506   Double_t cont = fHistAnode->GetCellContent(jc,ic);
1507   fPixArray->Add(new AliMUONPad (xc, yc, wx, wy, cont));
1508   used[(ic-1)*nx+jc-1] = kTRUE;
1509   AddBinSimple(fHistAnode, ic, jc);
1510   //fSplitter->AddBin(hist, ic, jc, 1, used, (TObjArray*)0); // recursive call
1511
1512   Int_t nPix = fPixArray->GetEntriesFast();
1513   Int_t npad = cluster.Multiplicity();
1514   
1515   for (Int_t i = 0; i < nPix; ++i) 
1516   {
1517     AliMUONPad* pixPtr = Pixel(i);
1518     pixPtr->SetSize(0,wx);
1519     pixPtr->SetSize(1,wy);
1520   }
1521
1522   // Pick up pads which overlap with found pixels
1523   for (Int_t i = 0; i < npad; ++i) 
1524   {
1525     //cluster.Pad(i)->SetStatus(-1);
1526     cluster.Pad(i)->SetStatus(fgkOver); // just the dirty trick
1527   }
1528   
1529   for (Int_t i = 0; i < nPix; ++i) 
1530   {
1531     AliMUONPad* pixPtr = Pixel(i);
1532     for (Int_t j = 0; j < npad; ++j) 
1533     {
1534       AliMUONPad* pad = cluster.Pad(j);
1535       //if (pad->Status() == 0) continue;
1536       if (pad->Status() == fgkZero) continue;
1537       if ( Overlap(*pad,*pixPtr) )
1538       {
1539         //pad->SetStatus(0);
1540         pad->SetStatus(fgkZero);
1541         if (fDebug) { cout << j << " "; pad->Print("full"); }
1542       }
1543     }
1544   }
1545
1546   delete [] used;
1547 }
1548
1549 //_____________________________________________________________________________
1550 void 
1551 AliMUONClusterFinderMLEM::AddBinSimple(TH2D *hist, Int_t ic, Int_t jc)
1552 {
1553   /// Add adjacent bins (+-1 in X and Y) to the cluster
1554   
1555   Int_t nx = hist->GetNbinsX();
1556   Int_t ny = hist->GetNbinsY();
1557   Double_t cont1, cont = hist->GetCellContent(jc,ic);
1558   AliMUONPad *pixPtr = 0;
1559   
1560   Int_t ie = TMath::Min(ic+1,ny), je = TMath::Min(jc+1,nx);
1561   for (Int_t i = TMath::Max(ic-1,1); i <= ie; ++i) {
1562     for (Int_t j = TMath::Max(jc-1,1); j <= je; ++j) {
1563       cont1 = hist->GetCellContent(j,i);
1564       if (cont1 > cont) continue;
1565       if (cont1 < fLowestPixelCharge) continue;
1566       pixPtr = new AliMUONPad (hist->GetXaxis()->GetBinCenter(j), 
1567                                hist->GetYaxis()->GetBinCenter(i), 0, 0, cont1);
1568       fPixArray->Add(pixPtr);
1569     }
1570   }
1571 }
1572
1573 //_____________________________________________________________________________
1574 AliMUONClusterFinderMLEM&  
1575 AliMUONClusterFinderMLEM::operator=(const AliMUONClusterFinderMLEM& rhs)
1576 {
1577 /// Protected assignement operator
1578
1579   if (this == &rhs) return *this;
1580
1581   AliFatal("Not implemented.");
1582     
1583   return *this;  
1584 }    
1585
1586 //_____________________________________________________________________________
1587 void AliMUONClusterFinderMLEM::AddVirtualPad(AliMUONCluster& cluster)
1588 {
1589   /// Add virtual pad (with small charge) to improve fit for clusters
1590   /// with number of pads == 2 per direction
1591   
1592   // Find out non-bending and bending planes
1593   Int_t nonb[2] = {1, 0}; // non-bending and bending cathodes
1594
1595   TVector2 dim0 = cluster.MinPadDimensions(0, 0, kTRUE);
1596   TVector2 dim1 = cluster.MinPadDimensions(1, 0, kTRUE);
1597   if (dim0.X() < dim1.X() - fgkDistancePrecision) {
1598     nonb[0] = 0;
1599     nonb[1] = 1;
1600   } 
1601
1602   Bool_t same = kFALSE;
1603   if (TMath::Abs(dim0.Y()-dim1.Y()) < fgkDistancePrecision) same = kTRUE; // the same pad size on both planes 
1604
1605   Long_t cn;
1606   Bool_t check[2] = {kFALSE, kFALSE};
1607   Int_t nxy[2];
1608   nxy[0] = nxy[1] = 0;
1609   for (Int_t inb = 0; inb < 2; ++inb) {
1610     cn = cluster.NofPads(nonb[inb], 0, kTRUE);
1611     if (inb == 0 && AliMp::PairFirst(cn) == 2) check[inb] = kTRUE; // check non-bending plane
1612     else if (inb == 1 && AliMp::PairSecond(cn) == 2) check[inb] = kTRUE; // check bending plane
1613     if (same) {
1614       nxy[0] = TMath::Max (nxy[0], AliMp::PairFirst(cn));
1615       nxy[1] = TMath::Max (nxy[1], AliMp::PairSecond(cn));
1616       if (inb == 0 && nxy[0] < 2) nonb[inb] = !nonb[inb];
1617       else if (inb == 1 && AliMp::PairSecond(cn) < 2) nonb[inb] = !nonb[inb];
1618     }
1619   }
1620   if (same) {
1621     if (nxy[0] > 2) check[0] = kFALSE;
1622     if (nxy[1] > 2) check[1] = kFALSE;
1623   }
1624   if (!check[0] && !check[1]) return;
1625
1626   for (Int_t inb = 0; inb < 2; ++inb) {
1627     if (!check[inb]) continue;
1628
1629     // Find pads with maximum and next to maximum charges 
1630     Int_t maxPads[2] = {-1, -1};
1631     Double_t amax[2] = {0};
1632     Int_t mult = cluster.Multiplicity();
1633     for (Int_t j = 0; j < mult; ++j) {
1634       AliMUONPad *pad = cluster.Pad(j);
1635       if (pad->Cathode() != nonb[inb]) continue;
1636       for (Int_t j2 = 0; j2 < 2; ++j2) {
1637         if (pad->Charge() > amax[j2]) {
1638           if (j2 == 0) { amax[1] = amax[0]; maxPads[1] = maxPads[0]; }
1639           amax[j2] = pad->Charge();
1640           maxPads[j2] = j;
1641           break;
1642         }
1643       }
1644     }
1645
1646     // Find min and max dimensions of the cluster
1647     Double_t limits[2] = {9999, -9999};
1648     for (Int_t j = 0; j < mult; ++j) {
1649       AliMUONPad *pad = cluster.Pad(j);
1650       if (pad->Cathode() != nonb[inb]) continue;
1651       if (pad->Coord(inb) < limits[0]) limits[0] = pad->Coord(inb);
1652       if (pad->Coord(inb) > limits[1]) limits[1] = pad->Coord(inb);
1653     }
1654
1655     // Loop over max and next to max pads
1656     Bool_t add = kFALSE;
1657     Int_t idirAdd = 0;
1658     for (Int_t j = 0; j < 2; ++j) {
1659       if (j == 1) {
1660         if (maxPads[j] < 0) continue;
1661         if (!add) break; 
1662         if (amax[1] / amax[0] < 0.5) break;
1663       }
1664       // Check if pad at the cluster limit
1665       AliMUONPad *pad = cluster.Pad(maxPads[j]);
1666       Int_t idir = 0;
1667       if (TMath::Abs(pad->Coord(inb)-limits[0]) < fgkDistancePrecision) idir = -1;
1668       else if (TMath::Abs(pad->Coord(inb)-limits[1]) < fgkDistancePrecision) idir = 1;
1669       else {
1670         //cout << " *** Pad not at the cluster limit: " << j << endl;
1671         break;
1672       }
1673       if (j == 1 && idir == idirAdd) break; // this pad has been already added
1674
1675       // Add pad (if it exists)
1676       TVector2 pos;
1677       if (inb == 0) pos.Set (pad->X() + idir * (pad->DX()+fgkDistancePrecision), pad->Y());
1678       else pos.Set (pad->X(), pad->Y() + idir * (pad->DY()+fgkDistancePrecision));
1679       //AliMpPad mppad = fkSegmentation[nonb[inb]]->PadByPosition(pos,kTRUE);
1680       AliMpPad mppad = fkSegmentation[nonb[inb]]->PadByPosition(pos.X(), pos.Y(),kFALSE);
1681       if (!mppad.IsValid()) continue; // non-existing pad
1682       AliMUONPad muonPad(fDetElemId, nonb[inb], mppad.GetIx(), mppad.GetIy(), 
1683                          mppad.GetPositionX(), mppad.GetPositionY(), 
1684                          mppad.GetDimensionX(), mppad.GetDimensionY(), 0);
1685       if (inb == 0) muonPad.SetCharge(TMath::Min (amax[j]/100, fLowestPadCharge));
1686       //else muonPad.SetCharge(TMath::Min (amax[j]/15, fgkZeroSuppression));
1687       else muonPad.SetCharge(TMath::Min (amax[j]/15, fLowestPadCharge));
1688       if (muonPad.Charge() < 2.0*fLowestPixelCharge) muonPad.SetCharge(2.0*fLowestPixelCharge);
1689       muonPad.SetReal(kFALSE);
1690       if (fDebug) printf(" ***** Add virtual pad in %d direction ***** %f %f %f %3d %3d %f %f \n",
1691                          inb, muonPad.Charge(), muonPad.X(), muonPad.Y(), muonPad.Ix(), 
1692                          muonPad.Iy(), muonPad.DX(), muonPad.DY());
1693       cluster.AddPad(muonPad); // add pad to the cluster
1694       add = kTRUE;
1695       idirAdd = idir;
1696     }
1697   }
1698 }
1699
1700 //_____________________________________________________________________________
1701 void AliMUONClusterFinderMLEM::PadsInXandY(AliMUONCluster& cluster,
1702                                            Int_t &nInX, Int_t &nInY) const
1703 {
1704   /// Find number of pads in X and Y-directions (excluding virtual ones and
1705   /// overflows)
1706
1707   //Int_t statusToTest = 1;
1708   Int_t statusToTest = fgkUseForFit;
1709   
1710   //if ( nInX < 0 ) statusToTest = 0;
1711   if ( nInX < 0 ) statusToTest = fgkZero;
1712        
1713   Bool_t mustMatch(kTRUE);
1714
1715   Long_t cn = cluster.NofPads(statusToTest,mustMatch);
1716   
1717   nInX = AliMp::PairFirst(cn);
1718   nInY = AliMp::PairSecond(cn);
1719 }
1720
1721 //_____________________________________________________________________________
1722 void AliMUONClusterFinderMLEM::RemovePixel(Int_t i)
1723 {
1724   /// Remove pixel at index i
1725   AliMUONPad* pixPtr = Pixel(i);
1726   fPixArray->RemoveAt(i); 
1727   delete pixPtr;
1728 }
1729
1730 //_____________________________________________________________________________
1731 void AliMUONClusterFinderMLEM::Simple(AliMUONCluster& cluster)
1732 {
1733 /// Process simple cluster (small number of pads) without EM-procedure
1734
1735   Int_t nForFit = 1, clustFit[1] = {0}, nfit;
1736   Double_t parOk[3] = {0.}; 
1737   TObjArray *clusters[1]; 
1738   clusters[0] = fPixArray;
1739
1740   AliDebug(1,Form("nPix=%d",fPixArray->GetLast()+1));
1741
1742   Int_t mult = cluster.Multiplicity();
1743   for (Int_t i = 0; i < mult; ++i) 
1744   {
1745     AliMUONPad* pad = cluster.Pad(i);
1746     /*
1747     if ( pad->IsSaturated()) 
1748     {
1749       pad->SetStatus(-9);
1750     }
1751     else 
1752     {
1753       pad->SetStatus(1);
1754     }
1755     */
1756     if (!pad->IsSaturated()) pad->SetStatus(fgkUseForFit);
1757   }
1758   nfit = fSplitter->Fit(cluster,1, nForFit, clustFit, clusters, parOk, fClusterList, fHistMlem);
1759 }
1760
1761 //_____________________________________________________________________________
1762 AliMUONPad* 
1763 AliMUONClusterFinderMLEM::Pixel(Int_t i) const
1764 {
1765   /// Returns pixel at index i
1766   return static_cast<AliMUONPad*>(fPixArray->UncheckedAt(i));
1767 }
1768
1769 //_____________________________________________________________________________
1770 void 
1771 AliMUONClusterFinderMLEM::Print(Option_t* what) const
1772 {
1773   /// printout
1774   TString swhat(what);
1775   swhat.ToLower();
1776   if ( swhat.Contains("precluster") )
1777   {
1778     if ( fPreCluster) fPreCluster->Print();
1779   }
1780 }
1781
1782 //_____________________________________________________________________________
1783 void 
1784 AliMUONClusterFinderMLEM::SetChargeHints(Double_t lowestPadCharge, Double_t lowestClusterCharge)
1785 {
1786   /// Set some thresholds we use later on in the algorithm
1787   fLowestPadCharge=lowestPadCharge;
1788   fLowestClusterCharge=lowestClusterCharge;
1789   fLowestPixelCharge=fLowestPadCharge/12.0; 
1790 }
1791
1792