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