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