]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONClusterFinderPeakFit.cxx
Coding violation fixes
[u/mrichter/AliRoot.git] / MUON / AliMUONClusterFinderPeakFit.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 AliMUONClusterFinderPeakFit
20 /// 
21 /// Clusterizer class based on simple peak finder
22 ///
23 /// Pre-clustering is handled by AliMUONPreClusterFinder
24 /// From a precluster a pixel array is built, and its local maxima are used
25 /// to get pads and make the fit with up to 3 hit candidates or compute pad 
26 /// centers of gravity for larger number of peaks.
27 ///
28 /// \author Laurent Aphecetche (for the "new" C++ structure) and 
29 /// Alexander Zinchenko, JINR Dubna, for the hardcore of it ;-)
30 //-----------------------------------------------------------------------------
31
32 #include "AliMUONClusterFinderPeakFit.h"
33 #include "AliMUONCluster.h"
34 #include "AliMUONConstants.h"
35 #include "AliMUONPad.h"
36 #include "AliMUONMathieson.h"
37
38 #include "AliMpDEManager.h"
39 #include "AliMpPad.h"
40 #include "AliMpVSegmentation.h"
41 #include "AliMpEncodePair.h"
42
43 #include "AliLog.h"
44 #include "AliRunLoader.h"
45 //#include "AliCodeTimer.h"
46
47 #include <Riostream.h>
48 #include <TH2.h>
49 #include <TVirtualFitter.h>
50 #include <TMath.h>
51 //#include <TCanvas.h>
52
53 /// \cond CLASSIMP
54 ClassImp(AliMUONClusterFinderPeakFit)
55 /// \endcond
56  
57 const Double_t AliMUONClusterFinderPeakFit::fgkZeroSuppression = 6; // average zero suppression value
58 //const Double_t AliMUONClusterFinderMLEM::fgkDistancePrecision = 1e-6; // (cm) used to check overlaps and so on
59 const Double_t AliMUONClusterFinderPeakFit::fgkDistancePrecision = 1e-3; // (cm) used to check overlaps and so on
60 const TVector2 AliMUONClusterFinderPeakFit::fgkIncreaseSize(-AliMUONClusterFinderPeakFit::fgkDistancePrecision,-AliMUONClusterFinderPeakFit::fgkDistancePrecision);
61 const TVector2 AliMUONClusterFinderPeakFit::fgkDecreaseSize(AliMUONClusterFinderPeakFit::fgkDistancePrecision,AliMUONClusterFinderPeakFit::fgkDistancePrecision);
62
63 // Status flags for pads
64 const Int_t AliMUONClusterFinderPeakFit::fgkZero = 0x0; ///< pad "basic" state
65 const Int_t AliMUONClusterFinderPeakFit::fgkMustKeep = 0x1; ///< do not kill (for pixels)
66 const Int_t AliMUONClusterFinderPeakFit::fgkUseForFit = 0x10; ///< should be used for fit
67 const Int_t AliMUONClusterFinderPeakFit::fgkOver = 0x100; ///< processing is over
68 const Int_t AliMUONClusterFinderPeakFit::fgkModified = 0x1000; ///< modified pad charge 
69 const Int_t AliMUONClusterFinderPeakFit::fgkCoupled = 0x10000; ///< coupled pad  
70
71 namespace
72 {
73   //_____________________________________________________________________________
74   Double_t Param2Coef(Int_t icand, Double_t coef, Double_t *par, Int_t nHits)
75   {
76     /// Extract hit contribution scale factor from fit parameters
77   
78     //Int_t nHits = TMath::Nint(par[8]);
79     if (nHits == 1) return 1.;
80     if (nHits == 2) return icand==0 ? par[2] : TMath::Max(1.-par[2],0.);
81     if (icand == 0) return par[2];
82     if (icand == 1) return TMath::Max((1.-par[2])*par[5], 0.);
83     return TMath::Max(1.-par[2]-coef,0.);
84   }
85
86   //___________________________________________________________________________
87   void 
88   FitFunction(Int_t& /*notused*/, Double_t* /*notused*/, 
89               Double_t& f, Double_t* par, 
90               Int_t /*notused*/)
91   {
92     /// Chi2 Function to minimize: Mathieson charge distribution in 2 dimensions
93     
94     TObjArray* userObjects = static_cast<TObjArray*>(TVirtualFitter::GetFitter()->GetObjectFit());
95     
96     AliMUONCluster* cluster = static_cast<AliMUONCluster*>(userObjects->At(0));
97     AliMUONMathieson* mathieson = static_cast<AliMUONMathieson*>(userObjects->At(1));
98     AliMUONClusterFinderPeakFit* finder = 
99       static_cast<AliMUONClusterFinderPeakFit*>(userObjects->At(2));
100     
101     f = 0.0;
102     Int_t nHits = finder->GetNMax(), npads = cluster->Multiplicity();
103     Double_t qTot = cluster->Charge(), coef = 0;
104     //if (cluster->Multiplicity(0) == 0 || cluster->Multiplicity(1) == 0) qTot *= 2.;
105     
106     for ( Int_t i = 0 ; i < npads; ++i )
107     {
108       AliMUONPad* pad = cluster->Pad(i);
109       // skip pads w/ saturation or other problem(s)
110       //if ( pad->Status() ) continue; 
111       if ( pad->IsSaturated() ) continue; 
112       Double_t charge = 0.;
113       for (Int_t j = 0; j < nHits; ++j) {
114         // Sum over hits
115         Int_t indx = 3 * j;
116         TVector2 lowerLeft = TVector2(par[indx],par[indx+1]) - pad->Position() - pad->Dimensions();
117         TVector2 upperRight(lowerLeft + pad->Dimensions()*2.0);
118         Double_t estimatedCharge = mathieson->IntXY(lowerLeft.X(),lowerLeft.Y(),
119                                                     upperRight.X(),upperRight.Y());
120         coef = Param2Coef(j, coef, par, nHits);
121         charge += estimatedCharge * coef;
122       }
123       charge *= qTot;
124       
125       Double_t delta = charge - pad->Charge();
126       delta *= delta;
127       delta /= pad->Charge();    
128       f += delta;
129     }
130     f /= (qTot/npads);
131     //cout << qTot << " " << par[0] << " " << par[1] << " " << f << endl;
132   }
133 }
134
135 //_____________________________________________________________________________
136 AliMUONClusterFinderPeakFit::AliMUONClusterFinderPeakFit(Bool_t plot, AliMUONVClusterFinder* clusterFinder)
137   : AliMUONVClusterFinder(),
138 fPreClusterFinder(clusterFinder),
139 fPreCluster(0x0),
140 fClusterList(),
141 fMathieson(0x0),
142 fEventNumber(0),
143 fDetElemId(-1),
144 fClusterNumber(0),
145 fNMax(0),
146 fHistAnode(0x0),
147 fPixArray(new TObjArray(20)),
148 fDebug(0),
149 fPlot(plot),
150 fNClusters(0),
151 fNAddVirtualPads(0)
152 {
153   /// Constructor
154  
155   fkSegmentation[1] = fkSegmentation[0] = 0x0; 
156
157   if (fPlot) fDebug = 1;
158 }
159
160 //_____________________________________________________________________________
161 AliMUONClusterFinderPeakFit::~AliMUONClusterFinderPeakFit()
162 {
163 /// Destructor
164   delete fPixArray; fPixArray = 0;
165   delete fPreClusterFinder;
166   AliInfo(Form("Total clusters %d AddVirtualPad needed %d",
167                fNClusters,fNAddVirtualPads));
168   delete fMathieson;
169 }
170
171 //_____________________________________________________________________________
172 Bool_t 
173 AliMUONClusterFinderPeakFit::Prepare(Int_t detElemId, TClonesArray* pads[2],
174                                      const AliMpArea& area, const AliMpVSegmentation* seg[2])
175 {
176   /// Prepare for clustering
177 //  AliCodeTimerAuto("",0)
178   
179   for ( Int_t i = 0; i < 2; ++i )
180   {
181     fkSegmentation[i] = seg[i];
182   }
183   
184   // Find out the DetElemId
185   fDetElemId = detElemId;
186   
187   // find out current event number, and reset the cluster number
188   AliRunLoader *runLoader = AliRunLoader::Instance();
189   fEventNumber = runLoader ? runLoader->GetEventNumber() : 0;
190   fClusterNumber = -1;
191   fClusterList.Delete();
192   
193   AliDebug(3,Form("EVT %d DE %d",fEventNumber,fDetElemId));
194   
195   AliMq::Station12Type stationType = AliMpDEManager::GetStation12Type(fDetElemId);
196   
197   Float_t kx3 = AliMUONConstants::SqrtKx3();
198   Float_t ky3 = AliMUONConstants::SqrtKy3();
199   Float_t pitch = AliMUONConstants::Pitch();
200   
201   if ( stationType == AliMq::kStation1 )
202   {
203     kx3 = AliMUONConstants::SqrtKx3St1();
204     ky3 = AliMUONConstants::SqrtKy3St1();
205     pitch = AliMUONConstants::PitchSt1();
206   }
207   
208   delete fMathieson;
209   fMathieson = new AliMUONMathieson;
210   
211   fMathieson->SetPitch(pitch);
212   fMathieson->SetSqrtKx3AndDeriveKx2Kx4(kx3);
213   fMathieson->SetSqrtKy3AndDeriveKy2Ky4(ky3);
214
215   if ( fPreClusterFinder->NeedSegmentation() )
216   {
217     return fPreClusterFinder->Prepare(detElemId,pads,area,seg);
218   }
219   else
220   {
221     return fPreClusterFinder->Prepare(detElemId,pads,area);
222   }
223 }
224
225 //_____________________________________________________________________________
226 AliMUONCluster* 
227 AliMUONClusterFinderPeakFit::NextCluster()
228 {
229   /// Return next cluster
230 //  AliCodeTimerAuto("",0)
231   
232   // if the list of clusters is not void, pick one from there
233   TObject* o = fClusterList.At(++fClusterNumber);
234   if ( o != 0x0 ) return static_cast<AliMUONCluster*>(o);
235   
236   //FIXME : at this point, must check whether we've used all the digits
237   //from precluster : if not, let the preclustering know about those unused
238   //digits, so it can reuse them
239   
240   // if the cluster list is exhausted, we need to go to the next
241   // pre-cluster and treat it
242
243   fClusterList.Delete(); // reset the list of clusters for this pre-cluster
244   fClusterNumber = -1; 
245     
246   fPreCluster = fPreClusterFinder->NextCluster();
247   
248   if (!fPreCluster)
249   {
250     // we are done
251     return 0x0;
252   }
253     
254   WorkOnPreCluster();
255
256   // WorkOnPreCluster may have used only part of the pads, so we check that
257   // now, and let the unused pads be reused by the preclustering...
258   
259   Int_t mult = fPreCluster->Multiplicity();
260   for ( Int_t i = 0; i < mult; ++i )
261   {
262     AliMUONPad* pad = fPreCluster->Pad(i);
263     if ( !pad->IsUsed() )
264     {
265       fPreClusterFinder->UsePad(*pad);
266     }
267   }
268   
269   return NextCluster();
270 }
271
272 //_____________________________________________________________________________
273 Bool_t
274 AliMUONClusterFinderPeakFit::WorkOnPreCluster()
275 {
276   /// Starting from a precluster, builds a pixel array, and then
277   /// extract clusters from this array
278   
279   //  AliCodeTimerAuto("",0)    
280
281   if (fDebug) {
282     cout << " *** Event # " << fEventNumber
283          << " det. elem.: " << fDetElemId << endl;
284     for (Int_t j = 0; j < fPreCluster->Multiplicity(); ++j) {
285       AliMUONPad* pad = fPreCluster->Pad(j);
286       printf(" bbb %3d %1d %8.4f %8.4f %8.4f %8.4f %6.1f %3d %3d %2d %1d %1d \n",
287              j, pad->Cathode(), pad->Coord(0), pad->Coord(1), pad->DX()*2, pad->DY()*2,
288              pad->Charge(), pad->Ix(), pad->Iy(), pad->Status(), pad->IsReal(), pad->IsSaturated());
289     }
290   }
291
292   AliMUONCluster* cluster = CheckPrecluster(*fPreCluster);
293   if (!cluster) return kFALSE;
294
295   BuildPixArray(*cluster);
296   
297   if ( fPixArray->GetLast() < 0 )
298   {
299     AliDebug(1,"No pixel for the above cluster");
300     delete cluster;
301     return kFALSE;
302   }
303   
304   Int_t nMax = 1, localMax[100], maxPos[100] = {0};
305   Double_t maxVal[100];
306   
307   nMax = FindLocalMaxima(fPixArray, localMax, maxVal); // find local maxima
308
309   if (nMax > 1) TMath::Sort(nMax, maxVal, maxPos, kTRUE); // in descending order
310   
311   if (nMax <= 3) {
312     FindClusterFit(*cluster, localMax, maxPos, nMax);
313   } else {
314     for (Int_t i = 0; i < nMax; ++i) 
315       {
316         FindClusterCOG(*cluster, localMax, maxPos[i]);
317       }
318   }
319
320   delete cluster;
321   if (fPlot == 0) {
322     delete fHistAnode;
323     fHistAnode = 0x0;
324   }
325   return kTRUE;
326 }
327
328 //_____________________________________________________________________________
329 Bool_t 
330 AliMUONClusterFinderPeakFit::Overlap(const AliMUONPad& pad, const AliMUONPad& pixel)
331 {
332   /// Check if the pad and the pixel overlaps
333
334   // make a fake pad from the pixel
335   AliMUONPad tmp(pad.DetElemId(),pad.Cathode(),pad.Ix(),pad.Iy(),
336                  pixel.Coord(0),pixel.Coord(1),
337                  pixel.Size(0),pixel.Size(1),0);
338   
339   return AliMUONPad::AreOverlapping(pad,tmp,fgkDecreaseSize);
340 }
341
342 //_____________________________________________________________________________
343 AliMUONCluster* 
344 AliMUONClusterFinderPeakFit::CheckPrecluster(const AliMUONCluster& origCluster)
345 {
346   /// Check precluster in order to attempt to simplify it (mostly for
347   /// two-cathode preclusters)
348     
349   //  AliCodeTimerAuto("",0)
350
351   // Disregard small clusters (leftovers from splitting or noise)
352   if ((origCluster.Multiplicity()==1 || origCluster.Multiplicity()==2) &&
353       origCluster.Charge(0)+origCluster.Charge(1) < 10) 
354   { 
355     return 0x0;
356   }
357
358   AliMUONCluster* cluster = new AliMUONCluster(origCluster);
359
360   AliDebug(2,"Start of CheckPreCluster=");
361   //StdoutToAliDebug(2,cluster->Print("full"));
362
363   AliMUONCluster* rv(0x0);
364   
365   if (cluster->Multiplicity(0) && cluster->Multiplicity(1))
366   { 
367     rv = CheckPreclusterTwoCathodes(cluster);
368   }
369   else
370   {
371     rv = cluster;
372   }
373   return rv;
374 }
375
376 //_____________________________________________________________________________
377 AliMUONCluster*
378 AliMUONClusterFinderPeakFit::CheckPreclusterTwoCathodes(AliMUONCluster* cluster)
379 {
380   /// Check two-cathode cluster
381   
382   Int_t npad = cluster->Multiplicity();
383   Int_t* flags = new Int_t[npad];
384   for (Int_t j = 0; j < npad; ++j) flags[j] = 0;
385   
386   // Check pad overlaps
387   for ( Int_t i = 0; i < npad; ++i) 
388   {
389     AliMUONPad* padi = cluster->Pad(i);
390     if ( padi->Cathode() != 0 ) continue;
391     for (Int_t j = i+1; j < npad; ++j) 
392     {
393       AliMUONPad* padj = cluster->Pad(j);
394       if ( padj->Cathode() != 1 ) continue;
395       if ( !AliMUONPad::AreOverlapping(*padi,*padj,fgkDecreaseSize) ) continue;
396       flags[i] = flags[j] = 1; // mark overlapped pads
397     } 
398   } 
399   
400   // Check if all pads overlap
401   Int_t nFlags=0;
402   for (Int_t i = 0; i < npad; ++i) 
403   {
404     if (!flags[i]) ++nFlags;
405   }
406   
407   if (nFlags > 0) 
408   {
409     // not all pads overlap.
410     if (fDebug) cout << " nFlags: " << nFlags << endl;
411     TObjArray toBeRemoved;
412     for (Int_t i = 0; i < npad; ++i) 
413     {
414       AliMUONPad* pad = cluster->Pad(i);
415       if (flags[i]) continue;
416       Int_t cath = pad->Cathode();
417       Int_t cath1 = TMath::Even(cath);
418       // Check for edge effect (missing pads on the _other_ cathode)
419       AliMpPad mpPad 
420         = fkSegmentation[cath1]->PadByPosition(pad->Position().X(),pad->Position().Y(),kFALSE);
421       if (!mpPad.IsValid()) continue;
422       //if (nFlags == 1 && pad->Charge() < fgkZeroSuppression * 3) continue;
423       if (nFlags == 1 && pad->Charge() < 20) continue;
424       AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
425                       fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),pad->Charge()));
426       toBeRemoved.AddLast(pad);
427       fPreCluster->Pad(i)->Release();
428     }
429     Int_t nRemove = toBeRemoved.GetEntriesFast();
430     for ( Int_t i = 0; i < nRemove; ++i )
431     {
432       cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
433     }
434   } 
435   
436   // Check correlations of cathode charges
437   if ( !cluster->IsSaturated() && cluster->ChargeAsymmetry() > 1 )
438   {
439     // big difference
440     Int_t cathode = cluster->MaxRawChargeCathode();
441     Int_t imin(-1);
442     Int_t imax(-1);
443     Double_t cmax(0);
444     Double_t cmin(1E9);
445     
446     // get min and max pad charges on the cathode opposite to the 
447     // max pad (given by MaxRawChargeCathode())
448     //
449     Int_t mult = cluster->Multiplicity();
450     for ( Int_t i = 0; i < mult; ++i )
451     {
452       AliMUONPad* pad = cluster->Pad(i);
453       if ( pad->Cathode() != cathode || !pad->IsReal() )
454       {
455         // only consider pads in the opposite cathode, and
456         // only consider real pads (i.e. exclude the virtual ones)
457         continue;
458       }
459       if ( pad->Charge() < cmin )
460       {
461         cmin = pad->Charge();
462         imin = i;
463         if (imax < 0) {
464           imax = imin;
465           cmax = cmin;
466         }
467       }
468       else if ( pad->Charge() > cmax )
469       {
470         cmax = pad->Charge();
471         imax = i;
472       }      
473     }
474     AliDebug(2,Form("Pad imin,imax %d,%d cmin,cmax %e,%e",
475                     imin,imax,cmin,cmax));
476     //
477     // arrange pads according to their distance to the max, normalized
478     // to the pad size
479     Double_t* dist = new Double_t[mult];
480     Double_t dxMin(1E9);
481     Double_t dyMin(1E9);
482     Double_t dmin(0);
483     
484     AliMUONPad* padmax = cluster->Pad(imax);
485     
486     for ( Int_t i = 0; i < mult; ++i )
487     {
488       dist[i] = 0.0;
489       if ( i == imax) continue;
490       AliMUONPad* pad = cluster->Pad(i);
491       if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
492       Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
493       Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
494       dist[i] = TMath::Sqrt(dx*dx+dy*dy);
495       if ( i == imin )
496       {
497         dmin = dist[i] + 1E-3; // distance to the pad with minimum charge
498         dxMin = dx;
499         dyMin = dy;
500       }      
501     }
502     
503     TMath::Sort(mult,dist,flags,kFALSE); // in ascending order
504     Double_t xmax(-1), distPrev(999);
505     TObjArray toBeRemoved;
506     
507     for ( Int_t i = 0; i < mult; ++i )
508     {
509       Int_t indx = flags[i];
510       AliMUONPad* pad = cluster->Pad(indx);
511       if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
512       if ( dist[indx] > dmin )
513       {
514         // farther than the minimum pad
515         Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
516         Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
517         dx *= dxMin;
518         dy *= dyMin;
519         if (dx >= 0 && dy >= 0) continue;
520         if (TMath::Abs(dx) > TMath::Abs(dy) && dx >= 0) continue;
521         if (TMath::Abs(dy) > TMath::Abs(dx) && dy >= 0) continue;        
522       }
523       if (dist[indx] > distPrev + 1) break; // overstepping empty pads
524       if ( pad->Charge() <= cmax || TMath::Abs(dist[indx]-xmax) < 1E-3 )
525       {
526         // release pad
527         if (TMath::Abs(dist[indx]-xmax) < 1.e-3) 
528         {
529           cmax = TMath::Max(pad->Charge(),cmax);
530         }
531         else
532         {
533           cmax = pad->Charge();
534         }
535         xmax = dist[indx];
536         distPrev = dist[indx];
537         AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
538                         fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),
539                         pad->Charge()));
540   
541         toBeRemoved.AddLast(pad);
542         fPreCluster->Pad(indx)->Release();
543       }
544     }
545     Int_t nRemove = toBeRemoved.GetEntriesFast();
546     for ( Int_t i = 0; i < nRemove; ++i )
547     {
548       cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
549     }    
550     delete[] dist;
551   } // if ( !cluster->IsSaturated() && 
552   
553   delete[] flags;
554   
555   AliDebug(2,"End of CheckPreClusterTwoCathodes=");
556   //StdoutToAliDebug(2,cluster->Print("full"));
557
558   return cluster;    
559 }
560
561 //_____________________________________________________________________________
562 void
563 AliMUONClusterFinderPeakFit::CheckOverlaps()
564 {
565   /// For debug only : check if some pixels overlap...
566   
567   Int_t nPix = fPixArray->GetLast()+1;
568   Int_t dummy(0);
569   
570   for ( Int_t i = 0; i < nPix; ++i )
571   {
572     AliMUONPad* pixelI = Pixel(i);
573     AliMUONPad pi(dummy,dummy,dummy,dummy,
574                   pixelI->Coord(0),pixelI->Coord(1),
575                   pixelI->Size(0),pixelI->Size(1),0.0);
576     
577     for ( Int_t j = i+1; j < nPix; ++j )
578     {
579       AliMUONPad* pixelJ = Pixel(j);
580       AliMUONPad pj(dummy,dummy,dummy,dummy,
581                     pixelJ->Coord(0),pixelJ->Coord(1),
582                     pixelJ->Size(0),pixelJ->Size(1),0.0);  
583       AliMpArea area;
584       
585       if ( AliMUONPad::AreOverlapping(pi,pj,fgkDecreaseSize,area) )
586       {
587         AliInfo(Form("The following 2 pixels (%d and %d) overlap !",i,j));
588         /*
589         StdoutToAliInfo(pixelI->Print();
590                         cout << " Surface = " << pixelI->Size(0)*pixelI->Size(1)*4 << endl;
591                         pixelJ->Print();
592                         cout << " Surface = " << pixelJ->Size(0)*pixelJ->Size(1)*4 << endl;
593                         cout << " Area surface = " << area.Dimensions().X()*area.Dimensions().Y()*4 << endl;
594                         cout << "-------" << endl;
595                         );
596         */        
597       }
598     }    
599   }
600 }
601
602 //_____________________________________________________________________________
603 void AliMUONClusterFinderPeakFit::BuildPixArray(AliMUONCluster& cluster)
604 {
605   /// Build pixel array 
606   
607   Int_t npad = cluster.Multiplicity();
608   if (npad<=0) 
609   {
610     AliWarning("Got no pad at all ?!");
611   }
612   
613   fPixArray->Delete();
614   BuildPixArrayOneCathode(cluster);
615   
616 //  StdoutToAliDebug(2,cout << "End of BuildPixelArray:" << endl;
617 //                   fPixArray->Print(););
618   //CheckOverlaps();//FIXME : this is for debug only. Remove it.
619 }
620
621 //_____________________________________________________________________________
622 void AliMUONClusterFinderPeakFit::BuildPixArrayOneCathode(AliMUONCluster& cluster)
623 {
624   /// Build the pixel array
625
626 //  AliDebug(2,Form("cluster.Multiplicity=%d",cluster.Multiplicity()));
627
628   TVector2 dim = cluster.MinPadDimensions (-1, kFALSE);
629   Double_t width[2] = {dim.X(), dim.Y()}, xy0[2];
630   Int_t found[2] = {0}, mult = cluster.Multiplicity();
631
632   for ( Int_t i = 0; i < mult; ++i) {
633     AliMUONPad* pad = cluster.Pad(i);
634     for (Int_t j = 0; j < 2; ++j) {
635       if (found[j] == 0 && TMath::Abs(pad->Size(j)-width[j]) < fgkDistancePrecision) { 
636         xy0[j] = pad->Coord(j);
637         found[j] = 1;
638       }
639     }
640     if (found[0] && found[1]) break;
641   }
642
643   Double_t min[2], max[2];
644   Int_t cath0 = 0, cath1 = 1;
645   if (cluster.Multiplicity(0) == 0) cath0 = 1;
646   else if (cluster.Multiplicity(1) == 0) cath1 = 0;
647
648   Double_t leftDownX, leftDownY;
649   cluster.Area(cath0).LeftDownCorner(leftDownX, leftDownY);
650   Double_t rightUpX, rightUpY; 
651   cluster.Area(cath0).RightUpCorner(rightUpX, rightUpY);
652   min[0] = leftDownX;
653   min[1] = leftDownY;
654   max[0] = rightUpX;
655   max[1] = rightUpY;
656   if (cath1 != cath0) {
657     cluster.Area(cath1).LeftDownCorner(leftDownX, leftDownY);
658     cluster.Area(cath1).RightUpCorner(rightUpX, rightUpY);
659     min[0] = TMath::Max (min[0], leftDownX);
660     min[1] = TMath::Max (min[1], leftDownY);
661     max[0] = TMath::Min (max[0], rightUpX);
662     max[1] = TMath::Min (max[1], rightUpY);
663   }
664
665   // Adjust limits
666   if (cath0 != cath1) {
667     TVector2 dim0 = cluster.MinPadDimensions (0, -1, kFALSE);
668     TVector2 dim1 = cluster.MinPadDimensions (1, -1, kFALSE);
669     if (TMath::Abs(dim0.Y()-dim1.Y()) < fgkDistancePrecision) {
670       // The same size of pads on both cathodes - check position
671       AliMUONPad* pad0 = cluster.Pad(0);
672       for ( Int_t i = 1; i < mult; ++i) {
673         AliMUONPad* pad = cluster.Pad(i);
674         if (pad->Cathode() == pad0->Cathode()) continue;
675         Double_t dist = TMath::Abs (pad0->Coord(1) - pad->Coord(1));
676         Double_t dd = dist - Int_t(dist/width[1]/2.) * width[1] * 2.;
677         if (TMath::Abs(dd/width[1]/2.-0.5) < fgkDistancePrecision) { 
678           // Half pad shift between cathodes
679           width[0] /= 2.; 
680           width[1] /= 2.;
681         }
682         break;
683       }
684     }
685   }
686
687   Int_t nbins[2];
688   for (Int_t i = 0; i < 2; ++i) {
689     Double_t dist = (min[i] - xy0[i]) / width[i] / 2;
690     if (TMath::Abs(dist) < 1.e-6) dist = -1.e-6;
691     min[i] = xy0[i] + (TMath::Nint(dist-TMath::Sign(1.e-6,dist)) 
692                        + TMath::Sign(0.5,dist)) * width[i] * 2;
693     nbins[i] = TMath::Nint ((max[i] - min[i]) / width[i] / 2);
694     if (nbins[i] == 0) ++nbins[i];
695     max[i] = min[i] + nbins[i] * width[i] * 2;
696     //cout << dist << " " << min[i] << " " << max[i] << " " << nbins[i] << endl;
697   }
698
699   // Book histogram
700   TH2D *hist1 = new TH2D ("Grid", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
701   TH2D *hist2 = new TH2D ("Entries", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
702   TAxis *xaxis = hist1->GetXaxis();
703   TAxis *yaxis = hist1->GetYaxis();
704
705   // Fill histogram
706   for ( Int_t i = 0; i < mult; ++i) {
707     AliMUONPad* pad = cluster.Pad(i);
708     Int_t ix0 = xaxis->FindBin(pad->X());
709     Int_t iy0 = yaxis->FindBin(pad->Y());
710     PadOverHist(0, ix0, iy0, pad, hist1, hist2);
711   }
712
713   // Store pixels
714   for (Int_t i = 1; i <= nbins[0]; ++i) {
715     Double_t x = xaxis->GetBinCenter(i);
716     for (Int_t j = 1; j <= nbins[1]; ++j) {
717       if (hist2->GetCellContent(i,j) < 0.1) continue;
718       if (cath0 != cath1) {
719         // Two-sided cluster
720         Double_t cont = hist2->GetCellContent(i,j);
721         if (cont < 999.) continue;
722         if (cont-Int_t(cont/1000.)*1000. < 0.5) continue;
723       }
724       //if (hist2->GetCellContent(i,j) < 1.1 && cluster.Multiplicity(0) && 
725       //  cluster.Multiplicity(1)) continue;
726       Double_t y = yaxis->GetBinCenter(j);
727       Double_t charge = hist1->GetCellContent(i,j);
728       AliMUONPad* pixPtr = new AliMUONPad(x, y, width[0], width[1], charge);
729       fPixArray->Add(pixPtr);
730     }  
731   }
732   /*
733   if (fPixArray->GetEntriesFast() == 1) {
734     // Split pixel into 2
735     AliMUONPad* pixPtr = static_cast<AliMUONPad*> (fPixArray->UncheckedAt(0));
736     pixPtr->SetSize(0,width[0]/2.);
737     pixPtr->Shift(0,-width[0]/4.);
738     pixPtr = new AliMUONPad(pixPtr->X()+width[0], pixPtr->Y(), width[0]/2., width[1], pixPtr->Charge());
739     fPixArray->Add(pixPtr);
740   }
741   */
742   //fPixArray->Print();
743   delete hist1;
744   delete hist2;
745 }
746
747 //_____________________________________________________________________________
748 void AliMUONClusterFinderPeakFit::PadOverHist(Int_t idir, Int_t ix0, Int_t iy0, AliMUONPad *pad,
749                                               TH2D *hist1, TH2D *hist2)
750 {
751   /// "Span" pad over histogram in the direction idir
752
753   TAxis *axis = idir == 0 ? hist1->GetXaxis() : hist1->GetYaxis();
754   Int_t nbins = axis->GetNbins(), cath = pad->Cathode();
755   Double_t bin = axis->GetBinWidth(1), amask = TMath::Power(1000.,cath*1.);
756
757   Int_t nbinPad = (Int_t)(pad->Size(idir)/bin*2+fgkDistancePrecision) + 1; // number of bins covered by pad
758
759   for (Int_t i = 0; i < nbinPad; ++i) {
760     Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
761     if (ixy > nbins) break;
762     Double_t lowEdge = axis->GetBinLowEdge(ixy);
763     if (lowEdge + fgkDistancePrecision > pad->Coord(idir) + pad->Size(idir)) break;
764     if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
765     else {
766       // Fill histogram
767       Double_t cont = pad->Charge();
768       if (hist2->GetCellContent(ix0, ixy) > 0.1) 
769         cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont) 
770              + TMath::Min (TMath::Max(hist1->GetCellContent(ix0, ixy),cont)*0.1, 10.);
771       hist1->SetCellContent(ix0, ixy, cont);
772       hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
773     }
774   }
775
776   for (Int_t i = -1; i > -nbinPad; --i) {
777     Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
778     if (ixy < 1) break;
779     Double_t upEdge = axis->GetBinUpEdge(ixy);
780     if (upEdge - fgkDistancePrecision < pad->Coord(idir) - pad->Size(idir)) break;
781     if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
782     else {
783       // Fill histogram
784       Double_t cont = pad->Charge();
785       if (hist2->GetCellContent(ix0, ixy) > 0.1) 
786         cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont)
787              + TMath::Min (TMath::Max(hist1->GetCellContent(ix0, ixy),cont)*0.1,10.);
788       hist1->SetCellContent(ix0, ixy, cont);
789       hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
790     }
791   }
792 }
793
794 //_____________________________________________________________________________
795 Int_t AliMUONClusterFinderPeakFit::FindLocalMaxima(TObjArray *pixArray, Int_t *localMax, Double_t *maxVal)
796 {
797 /// Find local maxima in pixel space 
798
799   AliDebug(1,Form("nPix=%d",pixArray->GetLast()+1));
800
801   //TH2D *hist = NULL;
802   //delete ((TH2D*) gROOT->FindObject("anode"));
803   //if (pixArray == fPixArray) hist = (TH2D*) gROOT->FindObject("anode");
804   //else { hist = (TH2D*) gROOT->FindObject("anode1"); cout << hist << endl; }
805   //if (hist) hist->Delete();
806   delete fHistAnode;
807  
808   Double_t xylim[4] = {999, 999, 999, 999};
809
810   Int_t nPix = pixArray->GetEntriesFast();
811   AliMUONPad *pixPtr = 0;
812   for (Int_t ipix = 0; ipix < nPix; ++ipix) {
813     pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
814     for (Int_t i = 0; i < 4; ++i) 
815          xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
816   }
817   for (Int_t i = 0; i < 4; ++i) xylim[i] -= pixPtr->Size(i/2); 
818
819   Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
820   Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
821   if (pixArray == fPixArray) fHistAnode = new TH2D("anode","anode",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
822   else fHistAnode = new TH2D("anode1","anode1",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
823   for (Int_t ipix = 0; ipix < nPix; ++ipix) {
824     pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
825     fHistAnode->Fill(pixPtr->Coord(0), pixPtr->Coord(1), pixPtr->Charge());
826   }
827 //  if (fDraw && pixArray == fPixArray) fDraw->DrawHist("c2", hist);
828
829   Int_t nMax = 0, indx, nxy = ny * nx;
830   Int_t *isLocalMax = new Int_t[nxy];
831   for (Int_t i = 0; i < nxy; ++i) isLocalMax[i] = 0; 
832
833   for (Int_t i = 1; i <= ny; ++i) {
834     indx = (i-1) * nx;
835     for (Int_t j = 1; j <= nx; ++j) {
836       if (fHistAnode->GetCellContent(j,i) < 0.5) continue;
837       //if (isLocalMax[indx+j-1] < 0) continue;
838       if (isLocalMax[indx+j-1] != 0) continue;
839       FlagLocalMax(fHistAnode, i, j, isLocalMax);
840     }
841   }
842
843   for (Int_t i = 1; i <= ny; ++i) {
844     indx = (i-1) * nx;
845     for (Int_t j = 1; j <= nx; ++j) {
846       if (isLocalMax[indx+j-1] > 0) { 
847         localMax[nMax] = indx + j - 1; 
848         maxVal[nMax++] = fHistAnode->GetCellContent(j,i);
849         if (nMax > 99) break;
850       }
851     }
852     if (nMax > 99) {
853       AliError(" Too many local maxima !!!");
854       break;
855     }
856   }
857   if (fDebug) cout << " Local max: " << nMax << endl;
858   delete [] isLocalMax; 
859   return nMax;
860 }
861
862 //_____________________________________________________________________________
863 void AliMUONClusterFinderPeakFit::FlagLocalMax(TH2D *hist, Int_t i, Int_t j, Int_t *isLocalMax)
864 {
865 /// Flag pixels (whether or not local maxima)
866
867   Int_t nx = hist->GetNbinsX();
868   Int_t ny = hist->GetNbinsY();
869   Int_t cont = TMath::Nint (hist->GetCellContent(j,i));
870   Int_t cont1 = 0, indx = (i-1)*nx+j-1, indx1 = 0, indx2 = 0;
871
872   Int_t ie = i + 2, je = j + 2;
873   for (Int_t i1 = i-1; i1 < ie; ++i1) {
874     if (i1 < 1 || i1 > ny) continue;
875     indx1 = (i1 - 1) * nx;
876     for (Int_t j1 = j-1; j1 < je; ++j1) {
877       if (j1 < 1 || j1 > nx) continue;
878       if (i == i1 && j == j1) continue;
879       indx2 = indx1 + j1 - 1;
880       cont1 = TMath::Nint (hist->GetCellContent(j1,i1));
881       if (cont < cont1) { isLocalMax[indx] = -1; return; }
882       else if (cont > cont1) isLocalMax[indx2] = -1;
883       else { // the same charge
884         isLocalMax[indx] = 1; 
885         if (isLocalMax[indx2] == 0) {
886           FlagLocalMax(hist, i1, j1, isLocalMax);
887           if (isLocalMax[indx2] < 0) { isLocalMax[indx] = -1; return; }
888           else isLocalMax[indx2] = -1;
889         }
890       } 
891     }
892   }
893   isLocalMax[indx] = 1; // local maximum
894 }
895
896 //_____________________________________________________________________________
897 void AliMUONClusterFinderPeakFit::FindClusterFit(AliMUONCluster& cluster, const Int_t *localMax, 
898                                                  const Int_t *maxPos, Int_t nMax)
899 {
900 /// Fit pad charge distribution with nMax hit hypothesis
901
902   //if (cluster.Multiplicity(0) == 0 || cluster.Multiplicity(1) == 0) cout << cluster.Multiplicity(0) << " " << cluster.Multiplicity(1) << " " << cluster.Charge() << " " << cluster.Charge(0) << " " << cluster.Charge(1) << " " << endl;
903
904   //TH2D *hist = (TH2D*) gROOT->FindObject("anode");
905   Int_t nx = fHistAnode->GetNbinsX();
906   //Int_t ny = hist->GetNbinsY();
907   Double_t xmin = fHistAnode->GetXaxis()->GetXmin(); //- hist->GetXaxis()->GetBinWidth(1);
908   Double_t xmax = fHistAnode->GetXaxis()->GetXmax(); //+ hist->GetXaxis()->GetBinWidth(1);
909   Double_t ymin = fHistAnode->GetYaxis()->GetXmin(); //- hist->GetYaxis()->GetBinWidth(1);
910   Double_t ymax = fHistAnode->GetYaxis()->GetXmax(); //+ hist->GetYaxis()->GetBinWidth(1);
911
912   TVirtualFitter* fitter = TVirtualFitter::Fitter(0,nMax*3);
913   fitter->Clear("");
914   fitter->SetFCN(FitFunction);
915
916   Float_t stepX = 0.01; // cm
917   Float_t stepY = 0.01; // cm
918   Float_t stepQ = 0.01; // 
919   
920   Double_t args[10] = {-1.}; // disable printout
921   
922   fitter->ExecuteCommand("SET PRINT",args,1);
923   fitter->ExecuteCommand("SET NOW",args,0); // no warnings
924   
925   Int_t indx = 0;
926   fNMax = nMax;
927   for (Int_t i = 0; i < nMax; ++i) {
928     Int_t ic = localMax[maxPos[i]] / nx + 1;
929     Int_t jc = localMax[maxPos[i]] % nx + 1;
930     Double_t yc = fHistAnode->GetYaxis()->GetBinCenter(ic);
931     Double_t xc = fHistAnode->GetXaxis()->GetBinCenter(jc);
932     indx = 3 * i;
933     fitter->SetParameter(indx,"Hit X position",xc,stepX,xmin,xmax);
934     fitter->SetParameter(indx+1,"Hit Y position",yc,stepY,ymin,ymax);
935     fitter->SetParameter(indx+2,"Hit contribution",0.6,stepQ,0.,1.);
936   }
937   fitter->SetParameter(indx+2,"Hit contribution",0.,0.,0,0);
938   //fitter->SetParameter(8,"Number of hits",nMax,0.,0,0);
939
940   TObjArray userObjects;
941   
942   userObjects.Add(&cluster);
943   userObjects.Add(fMathieson);
944   userObjects.Add(this);
945   
946   fitter->SetObjectFit(&userObjects);
947   
948   args[0] = 500.;
949   args[1] = 1.;
950   /*Int_t stat =*/ fitter->ExecuteCommand("MIGRAD",args,2);
951   //if (stat) { cout << " stat = " << stat << " " << fDetElemId << endl; /*exit(0);*/ }
952   //Int_t nvpar, nparx;
953   //Double_t amin, edm, errdef; 
954   //fitter->GetStats(amin, edm, errdef, nvpar, nparx);
955   //cout << amin << endl;
956
957   Double_t qTot = cluster.Charge(), par[9] = {0.}, err[9] = {0.}, coef = 0.;
958   //par[8] = nMax;
959   for (Int_t j = 0; j < nMax; ++j) {
960     indx = 3 * j;
961     par[indx+2] = fitter->GetParameter(indx+2);
962     coef = Param2Coef(j, coef, par, nMax);
963     par[indx] = fitter->GetParameter(indx);
964     par[indx+1] = fitter->GetParameter(indx+1);
965     err[indx] = fitter->GetParError(indx);
966     err[indx+1] = fitter->GetParError(indx+1);
967       
968     if ( coef*qTot >= 14 )
969     {
970       AliMUONCluster* cluster1 = new AliMUONCluster(cluster);
971       
972       cluster1->SetCharge(coef*qTot,coef*qTot);
973       
974       cluster1->SetPosition(TVector2(par[indx],par[indx+1]),TVector2(err[indx],err[indx+1]));
975       cluster1->SetChi2(0.);
976       
977       // FIXME: we miss some information in this cluster, as compared to 
978       // the original AddRawCluster code.
979       
980       AliDebug(2,Form("Adding RawCluster detElemId %4d mult %2d charge %5d (xl,yl)=(%9.6g,%9.6g)",
981                       fDetElemId,cluster1->Multiplicity(),(Int_t)cluster1->Charge(),
982                       cluster1->Position().X(),cluster1->Position().Y()));
983         
984       fClusterList.Add(cluster1);
985     }
986   }
987 }
988
989 //_____________________________________________________________________________
990 void AliMUONClusterFinderPeakFit::FindClusterCOG(AliMUONCluster& cluster, 
991                                                  const Int_t *localMax, Int_t iMax)
992 {
993 /// Find COG of pad charge distribution around local maximum \a iMax 
994
995   //TH2D *hist = (TH2D*) gROOT->FindObject("anode");
996   /* Just for check
997   TCanvas* c = new TCanvas("Anode","Anode",800,600);
998   c->cd();
999   hist->Draw("lego1Fb"); // debug
1000   c->Update();
1001   Int_t tmp;
1002   cin >> tmp;
1003   */
1004   Int_t nx = fHistAnode->GetNbinsX();
1005   //Int_t ny = hist->GetNbinsY();
1006   Int_t ic = localMax[iMax] / nx + 1;
1007   Int_t jc = localMax[iMax] % nx + 1;
1008
1009   // Get min pad dimensions for the precluster
1010   Int_t nSides = 2;
1011   if (cluster.Multiplicity(0) == 0 || cluster.Multiplicity(1) == 0) nSides = 1;
1012   TVector2 dim0 = cluster.MinPadDimensions(0, -1, kFALSE);
1013   TVector2 dim1 = cluster.MinPadDimensions(1, -1, kFALSE);
1014   //Double_t width[2][2] = {{dim0.X(), dim0.Y()},{dim1.X(),dim1.Y()}};
1015   Int_t nonb[2] = {1, 0}; // coordinate index vs cathode
1016   if (nSides == 1 || dim0.X() < dim1.X() - fgkDistancePrecision) {
1017     nonb[0] = 0;
1018     nonb[1] = 1;
1019   }
1020   Bool_t samex = kFALSE, samey = kFALSE;
1021   if (TMath::Abs(dim0.X()-dim1.X()) < fgkDistancePrecision) samex = kTRUE; // the same X pad size on both planes 
1022   if (TMath::Abs(dim0.Y()-dim1.Y()) < fgkDistancePrecision) samey = kTRUE; // the same Y pad size on both planes 
1023
1024   // Drop all pixels from the array - pick up only the ones from the cluster
1025   //fPixArray->Delete();
1026
1027   Double_t wx = fHistAnode->GetXaxis()->GetBinWidth(1)/2; 
1028   Double_t wy = fHistAnode->GetYaxis()->GetBinWidth(1)/2;  
1029   Double_t yc = fHistAnode->GetYaxis()->GetBinCenter(ic);
1030   Double_t xc = fHistAnode->GetXaxis()->GetBinCenter(jc);
1031   Double_t cont = fHistAnode->GetCellContent(jc,ic);
1032   AliMUONPad pixel(xc, yc, wx, wy, cont);
1033   if (fDebug) pixel.Print("full"); 
1034
1035   Int_t npad = cluster.Multiplicity();
1036   
1037   // Pick up pads which overlap with the maximum pixel and find pads with the max signal
1038   Double_t qMax[2] = {0}; 
1039   AliMUONPad *matrix[2][3] = {{0x0,0x0,0x0},{0x0,0x0,0x0}};
1040   for (Int_t j = 0; j < npad; ++j) 
1041   {
1042     AliMUONPad* pad = cluster.Pad(j);
1043     if ( Overlap(*pad,pixel) )
1044     {
1045       if (fDebug) { cout << j << " "; pad->Print("full"); }
1046       if (pad->Charge() > qMax[pad->Cathode()]) {
1047         qMax[pad->Cathode()] = pad->Charge();
1048         matrix[pad->Cathode()][1] = pad;
1049         if (nSides == 1) matrix[!pad->Cathode()][1] = pad;
1050       }
1051     }
1052   }
1053   //if (nSides == 2 && (matrix[0][1] == 0x0 || matrix[1][1] == 0x0)) return; // ???
1054
1055   // Find neighbours of maxima to have 3 pads per direction (if possible)
1056   for (Int_t j = 0; j < npad; ++j) 
1057   {
1058     AliMUONPad* pad = cluster.Pad(j);
1059     Int_t cath = pad->Cathode();
1060     if (pad == matrix[cath][1]) continue;
1061     Int_t nLoops = 3 - nSides;
1062
1063     for (Int_t k = 0; k < nLoops; ++k) {
1064       Int_t cath1 = cath;
1065       if (k) cath1 = !cath;
1066
1067       // Check the coordinate corresponding to the cathode (bending or non-bending case)
1068       Double_t dist = pad->Coord(nonb[cath1]) - matrix[cath][1]->Coord(nonb[cath1]);
1069       Double_t dir = TMath::Sign (1., dist);
1070       dist = TMath::Abs(dist) - pad->Size(nonb[cath1]) - matrix[cath][1]->Size(nonb[cath1]);
1071
1072       if (TMath::Abs(dist) < fgkDistancePrecision) {
1073         // Check the other coordinate
1074         dist = pad->Coord(!nonb[cath1]) - matrix[cath1][1]->Coord(!nonb[cath1]);
1075         if (TMath::Abs(dist) > 
1076             TMath::Max(pad->Size(!nonb[cath1]), matrix[cath1][1]->Size(!nonb[cath1])) - fgkDistancePrecision) break;
1077         Int_t idir = TMath::Nint (dir);
1078         if (matrix[cath1][1+idir] == 0x0) matrix[cath1][1+idir] = pad;
1079         else if (pad->Charge() > matrix[cath1][1+idir]->Charge()) matrix[cath1][1+idir] = pad; // diff. segmentation
1080         //cout << pad->Coord(nonb[cath1]) << " " << pad->Coord(!nonb[cath1]) << " " << pad->Size(nonb[cath1]) << " " << pad->Size(!nonb[cath1]) << " " << pad->Charge() << endl ;
1081         break;
1082       }
1083     }
1084   }
1085
1086   Double_t coord[2] = {0.}, qAver = 0.;
1087   for (Int_t i = 0; i < 2; ++i) {
1088     Double_t q = 0.;
1089     Double_t coordQ = 0.;
1090     Int_t cath = matrix[i][1]->Cathode();
1091     if (i && nSides == 1) cath = !cath;
1092     for (Int_t j = 0; j < 3; ++j) {
1093       if (matrix[i][j] == 0x0) continue;
1094       Double_t dq = matrix[i][j]->Charge();
1095       q += dq;
1096       coordQ += dq * matrix[i][j]->Coord(nonb[cath]);
1097       //coordQ += (matrix[i][j]->Charge() * matrix[i][j]->Coord(nonb[cath]));
1098     }
1099     coord[cath] = coordQ / q;
1100     qAver = TMath::Max (qAver, q);
1101   }
1102
1103   //qAver = TMath::Sqrt(qAver);
1104   if ( qAver >= 14 )
1105   {
1106     
1107     AliMUONCluster* cluster1 = new AliMUONCluster(cluster);
1108       
1109     cluster1->SetCharge(qAver,qAver);
1110     if (nonb[0] == 1) 
1111       cluster1->SetPosition(TVector2(coord[1],coord[0]),TVector2(0.,0.));
1112     else 
1113       cluster1->SetPosition(TVector2(coord[0],coord[1]),TVector2(0.,0.));
1114
1115     cluster1->SetChi2(0.);
1116       
1117     // FIXME: we miss some information in this cluster, as compared to 
1118     // the original AddRawCluster code.
1119       
1120     AliDebug(2,Form("Adding RawCluster detElemId %4d mult %2d charge %5d (xl,yl)=(%9.6g,%9.6g)",
1121                     fDetElemId,cluster1->Multiplicity(),(Int_t)cluster1->Charge(),
1122                     cluster1->Position().X(),cluster1->Position().Y()));
1123         
1124     fClusterList.Add(cluster1);
1125   }
1126 }
1127
1128 //_____________________________________________________________________________
1129 AliMUONClusterFinderPeakFit&  
1130 AliMUONClusterFinderPeakFit::operator=(const AliMUONClusterFinderPeakFit& rhs)
1131 {
1132 /// Protected assignement operator
1133
1134   if (this == &rhs) return *this;
1135
1136   AliFatal("Not implemented.");
1137     
1138   return *this;  
1139 }    
1140
1141 //_____________________________________________________________________________
1142 void AliMUONClusterFinderPeakFit::PadsInXandY(AliMUONCluster& cluster,
1143                                            Int_t &nInX, Int_t &nInY) const
1144 {
1145   /// Find number of pads in X and Y-directions (excluding virtual ones and
1146   /// overflows)
1147
1148   //Int_t statusToTest = 1;
1149   Int_t statusToTest = fgkUseForFit;
1150   
1151   //if ( nInX < 0 ) statusToTest = 0;
1152   if ( nInX < 0 ) statusToTest = fgkZero;
1153        
1154   Bool_t mustMatch(kTRUE);
1155
1156   Long_t cn = cluster.NofPads(statusToTest,mustMatch);
1157   
1158   nInX = AliMp::PairFirst(cn);
1159   nInY = AliMp::PairSecond(cn);
1160 }
1161
1162 //_____________________________________________________________________________
1163 void AliMUONClusterFinderPeakFit::RemovePixel(Int_t i)
1164 {
1165   /// Remove pixel at index i
1166   AliMUONPad* pixPtr = Pixel(i);
1167   fPixArray->RemoveAt(i); 
1168   delete pixPtr;
1169 }
1170
1171 //_____________________________________________________________________________
1172 AliMUONPad* 
1173 AliMUONClusterFinderPeakFit::Pixel(Int_t i) const
1174 {
1175   /// Returns pixel at index i
1176   return static_cast<AliMUONPad*>(fPixArray->UncheckedAt(i));
1177 }
1178
1179 //_____________________________________________________________________________
1180 void 
1181 AliMUONClusterFinderPeakFit::Print(Option_t* what) const
1182 {
1183   /// printout
1184   TString swhat(what);
1185   swhat.ToLower();
1186   if ( swhat.Contains("precluster") )
1187   {
1188     if ( fPreCluster) fPreCluster->Print();
1189   }
1190 }
1191
1192