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