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