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