]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONClusterFinderPeakFit.cxx
Raw2SDigits() restored and improved.
[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 "AliMpStationType.h"
41 #include "AliMpVSegmentation.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   fSegmentation[1] = fSegmentation[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("")
178   
179   for ( Int_t i = 0; i < 2; ++i )
180   {
181     fSegmentation[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   AliMp::StationType stationType = AliMpDEManager::GetStationType(fDetElemId);
196   
197   Float_t kx3 = AliMUONConstants::SqrtKx3();
198   Float_t ky3 = AliMUONConstants::SqrtKy3();
199   Float_t pitch = AliMUONConstants::Pitch();
200   
201   if ( stationType == AliMp::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("")
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   fPreCluster = fPreClusterFinder->NextCluster();
244   
245   if (!fPreCluster)
246   {
247     // we are done
248     return 0x0;
249   }
250     
251   fClusterList.Delete(); // reset the list of clusters for this pre-cluster
252   fClusterNumber = -1; 
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("")      
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("")
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 = fSegmentation[cath1]->PadByPosition(pad->Position(),kFALSE);
420       if (!mpPad.IsValid()) continue;
421       //if (nFlags == 1 && pad->Charge() < fgkZeroSuppression * 3) continue;
422       if (nFlags == 1 && pad->Charge() < 20) continue;
423       AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
424                       fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),pad->Charge()));
425       toBeRemoved.AddLast(pad);
426       fPreCluster->Pad(i)->Release();
427     }
428     Int_t nRemove = toBeRemoved.GetEntriesFast();
429     for ( Int_t i = 0; i < nRemove; ++i )
430     {
431       cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
432     }
433   } 
434   
435   // Check correlations of cathode charges
436   if ( !cluster->IsSaturated() && cluster->ChargeAsymmetry() > 1 )
437   {
438     // big difference
439     Int_t cathode = cluster->MaxRawChargeCathode();
440     Int_t imin(-1);
441     Int_t imax(-1);
442     Double_t cmax(0);
443     Double_t cmin(1E9);
444     
445     // get min and max pad charges on the cathode opposite to the 
446     // max pad (given by MaxRawChargeCathode())
447     //
448     Int_t mult = cluster->Multiplicity();
449     for ( Int_t i = 0; i < mult; ++i )
450     {
451       AliMUONPad* pad = cluster->Pad(i);
452       if ( pad->Cathode() != cathode || !pad->IsReal() )
453       {
454         // only consider pads in the opposite cathode, and
455         // only consider real pads (i.e. exclude the virtual ones)
456         continue;
457       }
458       if ( pad->Charge() < cmin )
459       {
460         cmin = pad->Charge();
461         imin = i;
462         if (imax < 0) {
463           imax = imin;
464           cmax = cmin;
465         }
466       }
467       else if ( pad->Charge() > cmax )
468       {
469         cmax = pad->Charge();
470         imax = i;
471       }      
472     }
473     AliDebug(2,Form("Pad imin,imax %d,%d cmin,cmax %e,%e",
474                     imin,imax,cmin,cmax));
475     //
476     // arrange pads according to their distance to the max, normalized
477     // to the pad size
478     Double_t* dist = new Double_t[mult];
479     Double_t dxMin(1E9);
480     Double_t dyMin(1E9);
481     Double_t dmin(0);
482     
483     AliMUONPad* padmax = cluster->Pad(imax);
484     
485     for ( Int_t i = 0; i < mult; ++i )
486     {
487       dist[i] = 0.0;
488       if ( i == imax) continue;
489       AliMUONPad* pad = cluster->Pad(i);
490       if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
491       Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
492       Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
493       dist[i] = TMath::Sqrt(dx*dx+dy*dy);
494       if ( i == imin )
495       {
496         dmin = dist[i] + 1E-3; // distance to the pad with minimum charge
497         dxMin = dx;
498         dyMin = dy;
499       }      
500     }
501     
502     TMath::Sort(mult,dist,flags,kFALSE); // in ascending order
503     Double_t xmax(-1), distPrev(999);
504     TObjArray toBeRemoved;
505     
506     for ( Int_t i = 0; i < mult; ++i )
507     {
508       Int_t indx = flags[i];
509       AliMUONPad* pad = cluster->Pad(indx);
510       if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
511       if ( dist[indx] > dmin )
512       {
513         // farther than the minimum pad
514         Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
515         Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
516         dx *= dxMin;
517         dy *= dyMin;
518         if (dx >= 0 && dy >= 0) continue;
519         if (TMath::Abs(dx) > TMath::Abs(dy) && dx >= 0) continue;
520         if (TMath::Abs(dy) > TMath::Abs(dx) && dy >= 0) continue;        
521       }
522       if (dist[indx] > distPrev + 1) break; // overstepping empty pads
523       if ( pad->Charge() <= cmax || TMath::Abs(dist[indx]-xmax) < 1E-3 )
524       {
525         // release pad
526         if (TMath::Abs(dist[indx]-xmax) < 1.e-3) 
527         {
528           cmax = TMath::Max(pad->Charge(),cmax);
529         }
530         else
531         {
532           cmax = pad->Charge();
533         }
534         xmax = dist[indx];
535         distPrev = dist[indx];
536         AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
537                         fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),
538                         pad->Charge()));
539   
540         toBeRemoved.AddLast(pad);
541         fPreCluster->Pad(indx)->Release();
542       }
543     }
544     Int_t nRemove = toBeRemoved.GetEntriesFast();
545     for ( Int_t i = 0; i < nRemove; ++i )
546     {
547       cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
548     }    
549     delete[] dist;
550   } // if ( !cluster->IsSaturated() && 
551   
552   delete[] flags;
553   
554   AliDebug(2,"End of CheckPreClusterTwoCathodes=");
555   //StdoutToAliDebug(2,cluster->Print("full"));
556
557   return cluster;    
558 }
559
560 //_____________________________________________________________________________
561 void
562 AliMUONClusterFinderPeakFit::CheckOverlaps()
563 {
564   /// For debug only : check if some pixels overlap...
565   
566   Int_t nPix = fPixArray->GetLast()+1;
567   Int_t dummy(0);
568   
569   for ( Int_t i = 0; i < nPix; ++i )
570   {
571     AliMUONPad* pixelI = Pixel(i);
572     AliMUONPad pi(dummy,dummy,dummy,dummy,
573                   pixelI->Coord(0),pixelI->Coord(1),
574                   pixelI->Size(0),pixelI->Size(1),0.0);
575     
576     for ( Int_t j = i+1; j < nPix; ++j )
577     {
578       AliMUONPad* pixelJ = Pixel(j);
579       AliMUONPad pj(dummy,dummy,dummy,dummy,
580                     pixelJ->Coord(0),pixelJ->Coord(1),
581                     pixelJ->Size(0),pixelJ->Size(1),0.0);  
582       AliMpArea area;
583       
584       if ( AliMUONPad::AreOverlapping(pi,pj,fgkDecreaseSize,area) )
585       {
586         AliInfo(Form("The following 2 pixels (%d and %d) overlap !",i,j));
587         /*
588         StdoutToAliInfo(pixelI->Print();
589                         cout << " Surface = " << pixelI->Size(0)*pixelI->Size(1)*4 << endl;
590                         pixelJ->Print();
591                         cout << " Surface = " << pixelJ->Size(0)*pixelJ->Size(1)*4 << endl;
592                         cout << " Area surface = " << area.Dimensions().X()*area.Dimensions().Y()*4 << endl;
593                         cout << "-------" << endl;
594                         );
595         */        
596       }
597     }    
598   }
599 }
600
601 //_____________________________________________________________________________
602 void AliMUONClusterFinderPeakFit::BuildPixArray(AliMUONCluster& cluster)
603 {
604   /// Build pixel array 
605   
606   Int_t npad = cluster.Multiplicity();
607   if (npad<=0) 
608   {
609     AliWarning("Got no pad at all ?!");
610   }
611   
612   fPixArray->Delete();
613   BuildPixArrayOneCathode(cluster);
614   
615 //  StdoutToAliDebug(2,cout << "End of BuildPixelArray:" << endl;
616 //                   fPixArray->Print(););
617   //CheckOverlaps();//FIXME : this is for debug only. Remove it.
618 }
619
620 //_____________________________________________________________________________
621 void AliMUONClusterFinderPeakFit::BuildPixArrayOneCathode(AliMUONCluster& cluster)
622 {
623   /// Build the pixel array
624
625 //  AliDebug(2,Form("cluster.Multiplicity=%d",cluster.Multiplicity()));
626
627   TVector2 dim = cluster.MinPadDimensions (-1, kFALSE);
628   Double_t width[2] = {dim.X(), dim.Y()}, xy0[2];
629   Int_t found[2] = {0}, mult = cluster.Multiplicity();
630
631   for ( Int_t i = 0; i < mult; ++i) {
632     AliMUONPad* pad = cluster.Pad(i);
633     for (Int_t j = 0; j < 2; ++j) {
634       if (found[j] == 0 && TMath::Abs(pad->Size(j)-width[j]) < fgkDistancePrecision) { 
635         xy0[j] = pad->Coord(j);
636         found[j] = 1;
637       }
638     }
639     if (found[0] && found[1]) break;
640   }
641
642   Double_t min[2], max[2];
643   Int_t cath0 = 0, cath1 = 1;
644   if (cluster.Multiplicity(0) == 0) cath0 = 1;
645   else if (cluster.Multiplicity(1) == 0) cath1 = 0;
646
647   TVector2 leftDown = cluster.Area(cath0).LeftDownCorner();
648   TVector2 rightUp = cluster.Area(cath0).RightUpCorner();
649   min[0] = leftDown.X();
650   min[1] = leftDown.Y();
651   max[0] = rightUp.X();
652   max[1] = rightUp.Y();
653   if (cath1 != cath0) {
654     leftDown = cluster.Area(cath1).LeftDownCorner();
655     rightUp = cluster.Area(cath1).RightUpCorner();
656     min[0] = TMath::Max (min[0], leftDown.X());
657     min[1] = TMath::Max (min[1], leftDown.Y());
658     max[0] = TMath::Min (max[0], rightUp.X());
659     max[1] = TMath::Min (max[1], rightUp.Y());
660   }
661
662   // Adjust limits
663   if (cath0 != cath1) {
664     TVector2 dim0 = cluster.MinPadDimensions (0, -1, kFALSE);
665     TVector2 dim1 = cluster.MinPadDimensions (1, -1, kFALSE);
666     if (TMath::Abs(dim0.Y()-dim1.Y()) < fgkDistancePrecision) {
667       // The same size of pads on both cathodes - check position
668       AliMUONPad* pad0 = cluster.Pad(0);
669       for ( Int_t i = 1; i < mult; ++i) {
670         AliMUONPad* pad = cluster.Pad(i);
671         if (pad->Cathode() == pad0->Cathode()) continue;
672         Double_t dist = TMath::Abs (pad0->Coord(1) - pad->Coord(1));
673         Double_t dd = dist - Int_t(dist/width[1]/2.) * width[1] * 2.;
674         if (TMath::Abs(dd/width[1]/2.-0.5) < fgkDistancePrecision) { 
675           // Half pad shift between cathodes
676           width[0] /= 2.; 
677           width[1] /= 2.;
678         }
679         break;
680       }
681     }
682   }
683
684   Int_t nbins[2];
685   for (Int_t i = 0; i < 2; ++i) {
686     Double_t dist = (min[i] - xy0[i]) / width[i] / 2;
687     if (TMath::Abs(dist) < 1.e-6) dist = -1.e-6;
688     min[i] = xy0[i] + (TMath::Nint(dist-TMath::Sign(1.e-6,dist)) 
689                        + TMath::Sign(0.5,dist)) * width[i] * 2;
690     nbins[i] = TMath::Nint ((max[i] - min[i]) / width[i] / 2);
691     if (nbins[i] == 0) ++nbins[i];
692     max[i] = min[i] + nbins[i] * width[i] * 2;
693     //cout << dist << " " << min[i] << " " << max[i] << " " << nbins[i] << endl;
694   }
695
696   // Book histogram
697   TH2D *hist1 = new TH2D ("Grid", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
698   TH2D *hist2 = new TH2D ("Entries", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
699   TAxis *xaxis = hist1->GetXaxis();
700   TAxis *yaxis = hist1->GetYaxis();
701
702   // Fill histogram
703   for ( Int_t i = 0; i < mult; ++i) {
704     AliMUONPad* pad = cluster.Pad(i);
705     Int_t ix0 = xaxis->FindBin(pad->X());
706     Int_t iy0 = yaxis->FindBin(pad->Y());
707     PadOverHist(0, ix0, iy0, pad, hist1, hist2);
708   }
709
710   // Store pixels
711   for (Int_t i = 1; i <= nbins[0]; ++i) {
712     Double_t x = xaxis->GetBinCenter(i);
713     for (Int_t j = 1; j <= nbins[1]; ++j) {
714       if (hist2->GetCellContent(i,j) < 0.1) continue;
715       if (cath0 != cath1) {
716         // Two-sided cluster
717         Double_t cont = hist2->GetCellContent(i,j);
718         if (cont < 999.) continue;
719         if (cont-Int_t(cont/1000.)*1000. < 0.5) continue;
720       }
721       //if (hist2->GetCellContent(i,j) < 1.1 && cluster.Multiplicity(0) && 
722       //  cluster.Multiplicity(1)) continue;
723       Double_t y = yaxis->GetBinCenter(j);
724       Double_t charge = hist1->GetCellContent(i,j);
725       AliMUONPad* pixPtr = new AliMUONPad(x, y, width[0], width[1], charge);
726       fPixArray->Add(pixPtr);
727     }  
728   }
729   /*
730   if (fPixArray->GetEntriesFast() == 1) {
731     // Split pixel into 2
732     AliMUONPad* pixPtr = static_cast<AliMUONPad*> (fPixArray->UncheckedAt(0));
733     pixPtr->SetSize(0,width[0]/2.);
734     pixPtr->Shift(0,-width[0]/4.);
735     pixPtr = new AliMUONPad(pixPtr->X()+width[0], pixPtr->Y(), width[0]/2., width[1], pixPtr->Charge());
736     fPixArray->Add(pixPtr);
737   }
738   */
739   //fPixArray->Print();
740   delete hist1;
741   delete hist2;
742 }
743
744 //_____________________________________________________________________________
745 void AliMUONClusterFinderPeakFit::PadOverHist(Int_t idir, Int_t ix0, Int_t iy0, AliMUONPad *pad,
746                                               TH2D *hist1, TH2D *hist2)
747 {
748   /// "Span" pad over histogram in the direction idir
749
750   TAxis *axis = idir == 0 ? hist1->GetXaxis() : hist1->GetYaxis();
751   Int_t nbins = axis->GetNbins(), cath = pad->Cathode();
752   Double_t bin = axis->GetBinWidth(1), amask = TMath::Power(1000.,cath*1.);
753
754   Int_t nbinPad = (Int_t)(pad->Size(idir)/bin*2+fgkDistancePrecision) + 1; // number of bins covered by pad
755
756   for (Int_t i = 0; i < nbinPad; ++i) {
757     Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
758     if (ixy > nbins) break;
759     Double_t lowEdge = axis->GetBinLowEdge(ixy);
760     if (lowEdge + fgkDistancePrecision > pad->Coord(idir) + pad->Size(idir)) break;
761     if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
762     else {
763       // Fill histogram
764       Double_t cont = pad->Charge();
765       if (hist2->GetCellContent(ix0, ixy) > 0.1) 
766         cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont) 
767              + TMath::Min (TMath::Max(hist1->GetCellContent(ix0, ixy),cont)*0.1, 10.);
768       hist1->SetCellContent(ix0, ixy, cont);
769       hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
770     }
771   }
772
773   for (Int_t i = -1; i > -nbinPad; --i) {
774     Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
775     if (ixy < 1) break;
776     Double_t upEdge = axis->GetBinUpEdge(ixy);
777     if (upEdge - fgkDistancePrecision < pad->Coord(idir) - pad->Size(idir)) break;
778     if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
779     else {
780       // Fill histogram
781       Double_t cont = pad->Charge();
782       if (hist2->GetCellContent(ix0, ixy) > 0.1) 
783         cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont)
784              + TMath::Min (TMath::Max(hist1->GetCellContent(ix0, ixy),cont)*0.1,10.);
785       hist1->SetCellContent(ix0, ixy, cont);
786       hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
787     }
788   }
789 }
790
791 //_____________________________________________________________________________
792 Int_t AliMUONClusterFinderPeakFit::FindLocalMaxima(TObjArray *pixArray, Int_t *localMax, Double_t *maxVal)
793 {
794 /// Find local maxima in pixel space 
795
796   AliDebug(1,Form("nPix=%d",pixArray->GetLast()+1));
797
798   //TH2D *hist = NULL;
799   //delete ((TH2D*) gROOT->FindObject("anode"));
800   //if (pixArray == fPixArray) hist = (TH2D*) gROOT->FindObject("anode");
801   //else { hist = (TH2D*) gROOT->FindObject("anode1"); cout << hist << endl; }
802   //if (hist) hist->Delete();
803   delete fHistAnode;
804  
805   Double_t xylim[4] = {999, 999, 999, 999};
806
807   Int_t nPix = pixArray->GetEntriesFast();
808   AliMUONPad *pixPtr = 0;
809   for (Int_t ipix = 0; ipix < nPix; ++ipix) {
810     pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
811     for (Int_t i = 0; i < 4; ++i) 
812          xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
813   }
814   for (Int_t i = 0; i < 4; ++i) xylim[i] -= pixPtr->Size(i/2); 
815
816   Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
817   Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
818   if (pixArray == fPixArray) fHistAnode = new TH2D("anode","anode",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
819   else fHistAnode = new TH2D("anode1","anode1",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
820   for (Int_t ipix = 0; ipix < nPix; ++ipix) {
821     pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
822     fHistAnode->Fill(pixPtr->Coord(0), pixPtr->Coord(1), pixPtr->Charge());
823   }
824 //  if (fDraw && pixArray == fPixArray) fDraw->DrawHist("c2", hist);
825
826   Int_t nMax = 0, indx, nxy = ny * nx;
827   Int_t *isLocalMax = new Int_t[nxy];
828   for (Int_t i = 0; i < nxy; ++i) isLocalMax[i] = 0; 
829
830   for (Int_t i = 1; i <= ny; ++i) {
831     indx = (i-1) * nx;
832     for (Int_t j = 1; j <= nx; ++j) {
833       if (fHistAnode->GetCellContent(j,i) < 0.5) continue;
834       //if (isLocalMax[indx+j-1] < 0) continue;
835       if (isLocalMax[indx+j-1] != 0) continue;
836       FlagLocalMax(fHistAnode, i, j, isLocalMax);
837     }
838   }
839
840   for (Int_t i = 1; i <= ny; ++i) {
841     indx = (i-1) * nx;
842     for (Int_t j = 1; j <= nx; ++j) {
843       if (isLocalMax[indx+j-1] > 0) { 
844         localMax[nMax] = indx + j - 1; 
845         maxVal[nMax++] = fHistAnode->GetCellContent(j,i);
846         if (nMax > 99) break;
847       }
848     }
849     if (nMax > 99) {
850       AliError(" Too many local maxima !!!");
851       break;
852     }
853   }
854   if (fDebug) cout << " Local max: " << nMax << endl;
855   delete [] isLocalMax; 
856   return nMax;
857 }
858
859 //_____________________________________________________________________________
860 void AliMUONClusterFinderPeakFit::FlagLocalMax(TH2D *hist, Int_t i, Int_t j, Int_t *isLocalMax)
861 {
862 /// Flag pixels (whether or not local maxima)
863
864   Int_t nx = hist->GetNbinsX();
865   Int_t ny = hist->GetNbinsY();
866   Int_t cont = TMath::Nint (hist->GetCellContent(j,i));
867   Int_t cont1 = 0, indx = (i-1)*nx+j-1, indx1 = 0, indx2 = 0;
868
869   Int_t ie = i + 2, je = j + 2;
870   for (Int_t i1 = i-1; i1 < ie; ++i1) {
871     if (i1 < 1 || i1 > ny) continue;
872     indx1 = (i1 - 1) * nx;
873     for (Int_t j1 = j-1; j1 < je; ++j1) {
874       if (j1 < 1 || j1 > nx) continue;
875       if (i == i1 && j == j1) continue;
876       indx2 = indx1 + j1 - 1;
877       cont1 = TMath::Nint (hist->GetCellContent(j1,i1));
878       if (cont < cont1) { isLocalMax[indx] = -1; return; }
879       else if (cont > cont1) isLocalMax[indx2] = -1;
880       else { // the same charge
881         isLocalMax[indx] = 1; 
882         if (isLocalMax[indx2] == 0) {
883           FlagLocalMax(hist, i1, j1, isLocalMax);
884           if (isLocalMax[indx2] < 0) { isLocalMax[indx] = -1; return; }
885           else isLocalMax[indx2] = -1;
886         }
887       } 
888     }
889   }
890   isLocalMax[indx] = 1; // local maximum
891 }
892
893 //_____________________________________________________________________________
894 void AliMUONClusterFinderPeakFit::FindClusterFit(AliMUONCluster& cluster, Int_t *localMax, 
895                                                  Int_t *maxPos, Int_t nMax)
896 {
897 /// Fit pad charge distribution with nMax hit hypothesis
898
899   //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;
900
901   //TH2D *hist = (TH2D*) gROOT->FindObject("anode");
902   Int_t nx = fHistAnode->GetNbinsX();
903   //Int_t ny = hist->GetNbinsY();
904   Double_t xmin = fHistAnode->GetXaxis()->GetXmin(); //- hist->GetXaxis()->GetBinWidth(1);
905   Double_t xmax = fHistAnode->GetXaxis()->GetXmax(); //+ hist->GetXaxis()->GetBinWidth(1);
906   Double_t ymin = fHistAnode->GetYaxis()->GetXmin(); //- hist->GetYaxis()->GetBinWidth(1);
907   Double_t ymax = fHistAnode->GetYaxis()->GetXmax(); //+ hist->GetYaxis()->GetBinWidth(1);
908
909   TVirtualFitter* fitter = TVirtualFitter::Fitter(0,nMax*3);
910   fitter->Clear("");
911   fitter->SetFCN(FitFunction);
912
913   Float_t stepX = 0.01; // cm
914   Float_t stepY = 0.01; // cm
915   Float_t stepQ = 0.01; // 
916   
917   Double_t args[10] = {-1.}; // disable printout
918   
919   fitter->ExecuteCommand("SET PRINT",args,1);
920   fitter->ExecuteCommand("SET NOW",args,0); // no warnings
921   
922   Int_t indx = 0;
923   fNMax = nMax;
924   for (Int_t i = 0; i < nMax; ++i) {
925     Int_t ic = localMax[maxPos[i]] / nx + 1;
926     Int_t jc = localMax[maxPos[i]] % nx + 1;
927     Double_t yc = fHistAnode->GetYaxis()->GetBinCenter(ic);
928     Double_t xc = fHistAnode->GetXaxis()->GetBinCenter(jc);
929     indx = 3 * i;
930     fitter->SetParameter(indx,"Hit X position",xc,stepX,xmin,xmax);
931     fitter->SetParameter(indx+1,"Hit Y position",yc,stepY,ymin,ymax);
932     fitter->SetParameter(indx+2,"Hit contribution",0.6,stepQ,0.,1.);
933   }
934   fitter->SetParameter(indx+2,"Hit contribution",0.,0.,0,0);
935   //fitter->SetParameter(8,"Number of hits",nMax,0.,0,0);
936
937   TObjArray userObjects;
938   
939   userObjects.Add(&cluster);
940   userObjects.Add(fMathieson);
941   userObjects.Add(this);
942   
943   fitter->SetObjectFit(&userObjects);
944   
945   args[0] = 500.;
946   args[1] = 1.;
947   /*Int_t stat =*/ fitter->ExecuteCommand("MIGRAD",args,2);
948   //if (stat) { cout << " stat = " << stat << " " << fDetElemId << endl; /*exit(0);*/ }
949   //Int_t nvpar, nparx;
950   //Double_t amin, edm, errdef; 
951   //fitter->GetStats(amin, edm, errdef, nvpar, nparx);
952   //cout << amin << endl;
953
954   Double_t qTot = cluster.Charge(), par[9] = {0.}, err[9] = {0.}, coef = 0.;
955   //par[8] = nMax;
956   for (Int_t j = 0; j < nMax; ++j) {
957     indx = 3 * j;
958     par[indx+2] = fitter->GetParameter(indx+2);
959     coef = Param2Coef(j, coef, par, nMax);
960     par[indx] = fitter->GetParameter(indx);
961     par[indx+1] = fitter->GetParameter(indx+1);
962     err[indx] = fitter->GetParError(indx);
963     err[indx+1] = fitter->GetParError(indx+1);
964       
965     if ( coef*qTot >= 14 )
966     {
967       AliMUONCluster* cluster1 = new AliMUONCluster(cluster);
968       
969       cluster1->SetCharge(coef*qTot,coef*qTot);
970       
971       cluster1->SetPosition(TVector2(par[indx],par[indx+1]),TVector2(err[indx],err[indx+1]));
972       cluster1->SetChi2(0.);
973       
974       // FIXME: we miss some information in this cluster, as compared to 
975       // the original AddRawCluster code.
976       
977       AliDebug(2,Form("Adding RawCluster detElemId %4d mult %2d charge %5d (xl,yl)=(%9.6g,%9.6g)",
978                       fDetElemId,cluster1->Multiplicity(),(Int_t)cluster1->Charge(),
979                       cluster1->Position().X(),cluster1->Position().Y()));
980         
981       fClusterList.Add(cluster1);
982     }
983   }
984 }
985
986 //_____________________________________________________________________________
987 void AliMUONClusterFinderPeakFit::FindClusterCOG(AliMUONCluster& cluster, 
988                                                  Int_t *localMax, Int_t iMax)
989 {
990 /// Find COG of pad charge distribution around local maximum \a iMax 
991
992   //TH2D *hist = (TH2D*) gROOT->FindObject("anode");
993   /* Just for check
994   TCanvas* c = new TCanvas("Anode","Anode",800,600);
995   c->cd();
996   hist->Draw("lego1Fb"); // debug
997   c->Update();
998   Int_t tmp;
999   cin >> tmp;
1000   */
1001   Int_t nx = fHistAnode->GetNbinsX();
1002   //Int_t ny = hist->GetNbinsY();
1003   Int_t ic = localMax[iMax] / nx + 1;
1004   Int_t jc = localMax[iMax] % nx + 1;
1005
1006   // Get min pad dimensions for the precluster
1007   Int_t nSides = 2;
1008   if (cluster.Multiplicity(0) == 0 || cluster.Multiplicity(1) == 0) nSides = 1;
1009   TVector2 dim0 = cluster.MinPadDimensions(0, -1, kFALSE);
1010   TVector2 dim1 = cluster.MinPadDimensions(1, -1, kFALSE);
1011   //Double_t width[2][2] = {{dim0.X(), dim0.Y()},{dim1.X(),dim1.Y()}};
1012   Int_t nonb[2] = {1, 0}; // coordinate index vs cathode
1013   if (nSides == 1 || dim0.X() < dim1.X() - fgkDistancePrecision) {
1014     nonb[0] = 0;
1015     nonb[1] = 1;
1016   }
1017   Bool_t samex = kFALSE, samey = kFALSE;
1018   if (TMath::Abs(dim0.X()-dim1.X()) < fgkDistancePrecision) samex = kTRUE; // the same X pad size on both planes 
1019   if (TMath::Abs(dim0.Y()-dim1.Y()) < fgkDistancePrecision) samey = kTRUE; // the same Y pad size on both planes 
1020
1021   // Drop all pixels from the array - pick up only the ones from the cluster
1022   //fPixArray->Delete();
1023
1024   Double_t wx = fHistAnode->GetXaxis()->GetBinWidth(1)/2; 
1025   Double_t wy = fHistAnode->GetYaxis()->GetBinWidth(1)/2;  
1026   Double_t yc = fHistAnode->GetYaxis()->GetBinCenter(ic);
1027   Double_t xc = fHistAnode->GetXaxis()->GetBinCenter(jc);
1028   Double_t cont = fHistAnode->GetCellContent(jc,ic);
1029   AliMUONPad* pixPtr = new AliMUONPad (xc, yc, wx, wy, cont);
1030   if (fDebug) pixPtr->Print("full"); 
1031
1032   Int_t npad = cluster.Multiplicity();
1033   
1034   // Pick up pads which overlap with the maximum pixel and find pads with the max signal
1035   Double_t qMax[2] = {0}; 
1036   AliMUONPad *matrix[2][3] = {{0x0,0x0,0x0},{0x0,0x0,0x0}};
1037   for (Int_t j = 0; j < npad; ++j) 
1038   {
1039     AliMUONPad* pad = cluster.Pad(j);
1040     if ( Overlap(*pad,*pixPtr) )
1041     {
1042       if (fDebug) { cout << j << " "; pad->Print("full"); }
1043       if (pad->Charge() > qMax[pad->Cathode()]) {
1044         qMax[pad->Cathode()] = pad->Charge();
1045         matrix[pad->Cathode()][1] = pad;
1046         if (nSides == 1) matrix[!pad->Cathode()][1] = pad;
1047       }
1048     }
1049   }
1050   //if (nSides == 2 && (matrix[0][1] == 0x0 || matrix[1][1] == 0x0)) return; // ???
1051
1052   // Find neighbours of maxima to have 3 pads per direction (if possible)
1053   for (Int_t j = 0; j < npad; ++j) 
1054   {
1055     AliMUONPad* pad = cluster.Pad(j);
1056     Int_t cath = pad->Cathode();
1057     if (pad == matrix[cath][1]) continue;
1058     Int_t nLoops = 3 - nSides;
1059
1060     for (Int_t k = 0; k < nLoops; ++k) {
1061       Int_t cath1 = cath;
1062       if (k) cath1 = !cath;
1063
1064       // Check the coordinate corresponding to the cathode (bending or non-bending case)
1065       Double_t dist = pad->Coord(nonb[cath1]) - matrix[cath][1]->Coord(nonb[cath1]);
1066       Double_t dir = TMath::Sign (1., dist);
1067       dist = TMath::Abs(dist) - pad->Size(nonb[cath1]) - matrix[cath][1]->Size(nonb[cath1]);
1068
1069       if (TMath::Abs(dist) < fgkDistancePrecision) {
1070         // Check the other coordinate
1071         dist = pad->Coord(!nonb[cath1]) - matrix[cath1][1]->Coord(!nonb[cath1]);
1072         if (TMath::Abs(dist) > 
1073             TMath::Max(pad->Size(!nonb[cath1]), matrix[cath1][1]->Size(!nonb[cath1])) - fgkDistancePrecision) break;
1074         Int_t idir = TMath::Nint (dir);
1075         if (matrix[cath1][1+idir] == 0x0) matrix[cath1][1+idir] = pad;
1076         else if (pad->Charge() > matrix[cath1][1+idir]->Charge()) matrix[cath1][1+idir] = pad; // diff. segmentation
1077         //cout << pad->Coord(nonb[cath1]) << " " << pad->Coord(!nonb[cath1]) << " " << pad->Size(nonb[cath1]) << " " << pad->Size(!nonb[cath1]) << " " << pad->Charge() << endl ;
1078         break;
1079       }
1080     }
1081   }
1082
1083   Double_t coord[2] = {0.}, qAver = 0.;
1084   for (Int_t i = 0; i < 2; ++i) {
1085     Double_t q = 0.;
1086     Double_t coordQ = 0.;
1087     Int_t cath = matrix[i][1]->Cathode();
1088     if (i && nSides == 1) cath = !cath;
1089     for (Int_t j = 0; j < 3; ++j) {
1090       if (matrix[i][j] == 0x0) continue;
1091       Double_t dq = matrix[i][j]->Charge();
1092       q += dq;
1093       coordQ += dq * matrix[i][j]->Coord(nonb[cath]);
1094       //coordQ += (matrix[i][j]->Charge() * matrix[i][j]->Coord(nonb[cath]));
1095     }
1096     coord[cath] = coordQ / q;
1097     qAver = TMath::Max (qAver, q);
1098   }
1099
1100   //qAver = TMath::Sqrt(qAver);
1101   if ( qAver >= 14 )
1102   {
1103     
1104     AliMUONCluster* cluster1 = new AliMUONCluster(cluster);
1105       
1106     cluster1->SetCharge(qAver,qAver);
1107     if (nonb[0] == 1) 
1108       cluster1->SetPosition(TVector2(coord[1],coord[0]),TVector2(0.,0.));
1109     else 
1110       cluster1->SetPosition(TVector2(coord[0],coord[1]),TVector2(0.,0.));
1111
1112     cluster1->SetChi2(0.);
1113       
1114     // FIXME: we miss some information in this cluster, as compared to 
1115     // the original AddRawCluster code.
1116       
1117     AliDebug(2,Form("Adding RawCluster detElemId %4d mult %2d charge %5d (xl,yl)=(%9.6g,%9.6g)",
1118                     fDetElemId,cluster1->Multiplicity(),(Int_t)cluster1->Charge(),
1119                     cluster1->Position().X(),cluster1->Position().Y()));
1120         
1121     fClusterList.Add(cluster1);
1122   }
1123 }
1124
1125 //_____________________________________________________________________________
1126 AliMUONClusterFinderPeakFit&  
1127 AliMUONClusterFinderPeakFit::operator=(const AliMUONClusterFinderPeakFit& rhs)
1128 {
1129 /// Protected assignement operator
1130
1131   if (this == &rhs) return *this;
1132
1133   AliFatal("Not implemented.");
1134     
1135   return *this;  
1136 }    
1137
1138 //_____________________________________________________________________________
1139 void AliMUONClusterFinderPeakFit::PadsInXandY(AliMUONCluster& cluster,
1140                                            Int_t &nInX, Int_t &nInY) const
1141 {
1142   /// Find number of pads in X and Y-directions (excluding virtual ones and
1143   /// overflows)
1144
1145   //Int_t statusToTest = 1;
1146   Int_t statusToTest = fgkUseForFit;
1147   
1148   //if ( nInX < 0 ) statusToTest = 0;
1149   if ( nInX < 0 ) statusToTest = fgkZero;
1150        
1151   Bool_t mustMatch(kTRUE);
1152
1153   AliMpIntPair cn = cluster.NofPads(statusToTest,mustMatch);
1154   
1155   nInX = cn.GetFirst();
1156   nInY = cn.GetSecond();
1157 }
1158
1159 //_____________________________________________________________________________
1160 void AliMUONClusterFinderPeakFit::RemovePixel(Int_t i)
1161 {
1162   /// Remove pixel at index i
1163   AliMUONPad* pixPtr = Pixel(i);
1164   fPixArray->RemoveAt(i); 
1165   delete pixPtr;
1166 }
1167
1168 //_____________________________________________________________________________
1169 AliMUONPad* 
1170 AliMUONClusterFinderPeakFit::Pixel(Int_t i) const
1171 {
1172   /// Returns pixel at index i
1173   return static_cast<AliMUONPad*>(fPixArray->UncheckedAt(i));
1174 }
1175
1176 //_____________________________________________________________________________
1177 void 
1178 AliMUONClusterFinderPeakFit::Print(Option_t* what) const
1179 {
1180   /// printout
1181   TString swhat(what);
1182   swhat.ToLower();
1183   if ( swhat.Contains("precluster") )
1184   {
1185     if ( fPreCluster) fPreCluster->Print();
1186   }
1187 }
1188
1189