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