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