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