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