1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 //-----------------------------------------------------------------------------
19 /// \class AliMUONClusterFinderPeakCOG
21 /// Clusterizer class based on simple peak finder
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 compute pad center of gravity.
27 /// \author Laurent Aphecetche (for the "new" C++ structure) and
28 /// Alexander Zinchenko, JINR Dubna, for the hardcore of it ;-)
29 //-----------------------------------------------------------------------------
31 #include "AliMUONClusterFinderPeakCOG.h"
32 #include "AliMUONCluster.h"
33 #include "AliMUONPad.h"
36 #include "AliMpVSegmentation.h"
39 #include "AliRunLoader.h"
40 //#include "AliCodeTimer.h"
42 #include <Riostream.h>
44 //#include <TCanvas.h>
48 ClassImp(AliMUONClusterFinderPeakCOG)
51 const Double_t AliMUONClusterFinderPeakCOG::fgkZeroSuppression = 6; // average zero suppression value
52 //const Double_t AliMUONClusterFinderMLEM::fgkDistancePrecision = 1e-6; // (cm) used to check overlaps and so on
53 const Double_t AliMUONClusterFinderPeakCOG::fgkDistancePrecision = 1e-3; // (cm) used to check overlaps and so on
54 const TVector2 AliMUONClusterFinderPeakCOG::fgkIncreaseSize(-AliMUONClusterFinderPeakCOG::fgkDistancePrecision,-AliMUONClusterFinderPeakCOG::fgkDistancePrecision);
55 const TVector2 AliMUONClusterFinderPeakCOG::fgkDecreaseSize(AliMUONClusterFinderPeakCOG::fgkDistancePrecision,AliMUONClusterFinderPeakCOG::fgkDistancePrecision);
57 // Status flags for pads
58 const Int_t AliMUONClusterFinderPeakCOG::fgkZero = 0x0; ///< pad "basic" state
59 const Int_t AliMUONClusterFinderPeakCOG::fgkMustKeep = 0x1; ///< do not kill (for pixels)
60 const Int_t AliMUONClusterFinderPeakCOG::fgkUseForFit = 0x10; ///< should be used for fit
61 const Int_t AliMUONClusterFinderPeakCOG::fgkOver = 0x100; ///< processing is over
62 const Int_t AliMUONClusterFinderPeakCOG::fgkModified = 0x1000; ///< modified pad charge
63 const Int_t AliMUONClusterFinderPeakCOG::fgkCoupled = 0x10000; ///< coupled pad
65 //_____________________________________________________________________________
66 AliMUONClusterFinderPeakCOG::AliMUONClusterFinderPeakCOG(Bool_t plot, AliMUONVClusterFinder* clusterFinder)
67 : AliMUONVClusterFinder(),
68 fPreClusterFinder(clusterFinder),
75 fPixArray(new TObjArray(20)),
83 fSegmentation[1] = fSegmentation[0] = 0x0;
85 if (fPlot) fDebug = 1;
88 //_____________________________________________________________________________
89 AliMUONClusterFinderPeakCOG::~AliMUONClusterFinderPeakCOG()
92 delete fPixArray; fPixArray = 0;
93 delete fPreClusterFinder;
94 AliInfo(Form("Total clusters %d AddVirtualPad needed %d",
95 fNClusters,fNAddVirtualPads));
98 //_____________________________________________________________________________
100 AliMUONClusterFinderPeakCOG::Prepare(Int_t detElemId, TClonesArray* pads[2],
101 const AliMpArea& area, const AliMpVSegmentation* seg[2])
103 /// Prepare for clustering
104 // AliCodeTimerAuto("")
106 for ( Int_t i = 0; i < 2; ++i )
108 fSegmentation[i] = seg[i];
111 // Find out the DetElemId
112 fDetElemId = detElemId;
114 // find out current event number, and reset the cluster number
115 AliRunLoader *runLoader = AliRunLoader::GetRunLoader();
116 fEventNumber = runLoader ? runLoader->GetEventNumber() : 0;
118 fClusterList.Delete();
120 AliDebug(3,Form("EVT %d DE %d",fEventNumber,fDetElemId));
122 if ( fPreClusterFinder->NeedSegmentation() )
124 return fPreClusterFinder->Prepare(detElemId,pads,area,seg);
128 return fPreClusterFinder->Prepare(detElemId,pads,area);
132 //_____________________________________________________________________________
134 AliMUONClusterFinderPeakCOG::NextCluster()
136 /// Return next cluster
137 // AliCodeTimerAuto("")
139 // if the list of clusters is not void, pick one from there
140 TObject* o = fClusterList.At(++fClusterNumber);
141 if ( o != 0x0 ) return static_cast<AliMUONCluster*>(o);
143 //FIXME : at this point, must check whether we've used all the digits
144 //from precluster : if not, let the preclustering know about those unused
145 //digits, so it can reuse them
147 // if the cluster list is exhausted, we need to go to the next
148 // pre-cluster and treat it
150 fPreCluster = fPreClusterFinder->NextCluster();
158 fClusterList.Delete(); // reset the list of clusters for this pre-cluster
163 // WorkOnPreCluster may have used only part of the pads, so we check that
164 // now, and let the unused pads be reused by the preclustering...
166 Int_t mult = fPreCluster->Multiplicity();
167 for ( Int_t i = 0; i < mult; ++i )
169 AliMUONPad* pad = fPreCluster->Pad(i);
170 if ( !pad->IsUsed() )
172 fPreClusterFinder->UsePad(*pad);
176 return NextCluster();
179 //_____________________________________________________________________________
181 AliMUONClusterFinderPeakCOG::WorkOnPreCluster()
183 /// Starting from a precluster, builds a pixel array, and then
184 /// extract clusters from this array
186 // AliCodeTimerAuto("")
189 cout << " *** Event # " << fEventNumber
190 << " det. elem.: " << fDetElemId << endl;
191 for (Int_t j = 0; j < fPreCluster->Multiplicity(); ++j) {
192 AliMUONPad* pad = fPreCluster->Pad(j);
193 printf(" bbb %3d %1d %8.4f %8.4f %8.4f %8.4f %6.1f %3d %3d %2d %1d %1d \n",
194 j, pad->Cathode(), pad->Coord(0), pad->Coord(1), pad->DX()*2, pad->DY()*2,
195 pad->Charge(), pad->Ix(), pad->Iy(), pad->Status(), pad->IsReal(), pad->IsSaturated());
199 AliMUONCluster* cluster = CheckPrecluster(*fPreCluster);
200 if (!cluster) return kFALSE;
202 BuildPixArray(*cluster);
204 if ( fPixArray->GetLast() < 0 )
206 AliDebug(1,"No pixel for the above cluster");
211 Int_t nMax = 1, localMax[100], maxPos[100] = {0};
212 Double_t maxVal[100];
214 nMax = FindLocalMaxima(fPixArray, localMax, maxVal); // find local maxima
216 if (nMax > 1) TMath::Sort(nMax, maxVal, maxPos, kTRUE); // in descending order
218 for (Int_t i = 0; i < nMax; ++i)
220 FindCluster(*cluster,localMax, maxPos[i]);
231 //_____________________________________________________________________________
233 AliMUONClusterFinderPeakCOG::Overlap(const AliMUONPad& pad, const AliMUONPad& pixel)
235 /// Check if the pad and the pixel overlaps
237 // make a fake pad from the pixel
238 AliMUONPad tmp(pad.DetElemId(),pad.Cathode(),pad.Ix(),pad.Iy(),
239 pixel.Coord(0),pixel.Coord(1),
240 pixel.Size(0),pixel.Size(1),0);
242 return AliMUONPad::AreOverlapping(pad,tmp,fgkDecreaseSize);
245 //_____________________________________________________________________________
247 AliMUONClusterFinderPeakCOG::CheckPrecluster(const AliMUONCluster& origCluster)
249 /// Check precluster in order to attempt to simplify it (mostly for
250 /// two-cathode preclusters)
252 // AliCodeTimerAuto("")
254 // Disregard small clusters (leftovers from splitting or noise)
255 if ((origCluster.Multiplicity()==1 || origCluster.Multiplicity()==2) &&
256 origCluster.Charge(0)+origCluster.Charge(1) < 10)
261 AliMUONCluster* cluster = new AliMUONCluster(origCluster);
263 AliDebug(2,"Start of CheckPreCluster=");
264 //StdoutToAliDebug(2,cluster->Print("full"));
266 AliMUONCluster* rv(0x0);
268 if (cluster->Multiplicity(0) && cluster->Multiplicity(1))
270 rv = CheckPreclusterTwoCathodes(cluster);
279 //_____________________________________________________________________________
281 AliMUONClusterFinderPeakCOG::CheckPreclusterTwoCathodes(AliMUONCluster* cluster)
283 /// Check two-cathode cluster
285 Int_t npad = cluster->Multiplicity();
286 Int_t* flags = new Int_t[npad];
287 for (Int_t j = 0; j < npad; ++j) flags[j] = 0;
289 // Check pad overlaps
290 for ( Int_t i = 0; i < npad; ++i)
292 AliMUONPad* padi = cluster->Pad(i);
293 if ( padi->Cathode() != 0 ) continue;
294 for (Int_t j = i+1; j < npad; ++j)
296 AliMUONPad* padj = cluster->Pad(j);
297 if ( padj->Cathode() != 1 ) continue;
298 if ( !AliMUONPad::AreOverlapping(*padi,*padj,fgkDecreaseSize) ) continue;
299 flags[i] = flags[j] = 1; // mark overlapped pads
303 // Check if all pads overlap
305 for (Int_t i = 0; i < npad; ++i)
307 if (!flags[i]) ++nFlags;
312 // not all pads overlap.
313 if (fDebug) cout << " nFlags: " << nFlags << endl;
314 TObjArray toBeRemoved;
315 for (Int_t i = 0; i < npad; ++i)
317 AliMUONPad* pad = cluster->Pad(i);
318 if (flags[i]) continue;
319 Int_t cath = pad->Cathode();
320 Int_t cath1 = TMath::Even(cath);
321 // Check for edge effect (missing pads on the _other_ cathode)
322 AliMpPad mpPad = fSegmentation[cath1]->PadByPosition(pad->Position(),kFALSE);
323 if (!mpPad.IsValid()) continue;
324 //if (nFlags == 1 && pad->Charge() < fgkZeroSuppression * 3) continue;
325 if (nFlags == 1 && pad->Charge() < 20) continue;
326 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
327 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),pad->Charge()));
328 toBeRemoved.AddLast(pad);
329 fPreCluster->Pad(i)->Release();
331 Int_t nRemove = toBeRemoved.GetEntriesFast();
332 for ( Int_t i = 0; i < nRemove; ++i )
334 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
338 // Check correlations of cathode charges
339 if ( !cluster->IsSaturated() && cluster->ChargeAsymmetry() > 1 )
342 Int_t cathode = cluster->MaxRawChargeCathode();
348 // get min and max pad charges on the cathode opposite to the
349 // max pad (given by MaxRawChargeCathode())
351 Int_t mult = cluster->Multiplicity();
352 for ( Int_t i = 0; i < mult; ++i )
354 AliMUONPad* pad = cluster->Pad(i);
355 if ( pad->Cathode() != cathode || !pad->IsReal() )
357 // only consider pads in the opposite cathode, and
358 // only consider real pads (i.e. exclude the virtual ones)
361 if ( pad->Charge() < cmin )
363 cmin = pad->Charge();
370 else if ( pad->Charge() > cmax )
372 cmax = pad->Charge();
376 AliDebug(2,Form("Pad imin,imax %d,%d cmin,cmax %e,%e",
377 imin,imax,cmin,cmax));
379 // arrange pads according to their distance to the max, normalized
381 Double_t* dist = new Double_t[mult];
386 AliMUONPad* padmax = cluster->Pad(imax);
388 for ( Int_t i = 0; i < mult; ++i )
391 if ( i == imax) continue;
392 AliMUONPad* pad = cluster->Pad(i);
393 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
394 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
395 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
396 dist[i] = TMath::Sqrt(dx*dx+dy*dy);
399 dmin = dist[i] + 1E-3; // distance to the pad with minimum charge
405 TMath::Sort(mult,dist,flags,kFALSE); // in ascending order
406 Double_t xmax(-1), distPrev(999);
407 TObjArray toBeRemoved;
409 for ( Int_t i = 0; i < mult; ++i )
411 Int_t indx = flags[i];
412 AliMUONPad* pad = cluster->Pad(indx);
413 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
414 if ( dist[indx] > dmin )
416 // farther than the minimum pad
417 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
418 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
421 if (dx >= 0 && dy >= 0) continue;
422 if (TMath::Abs(dx) > TMath::Abs(dy) && dx >= 0) continue;
423 if (TMath::Abs(dy) > TMath::Abs(dx) && dy >= 0) continue;
425 if (dist[indx] > distPrev + 1) break; // overstepping empty pads
426 if ( pad->Charge() <= cmax || TMath::Abs(dist[indx]-xmax) < 1E-3 )
429 if (TMath::Abs(dist[indx]-xmax) < 1.e-3)
431 cmax = TMath::Max(pad->Charge(),cmax);
435 cmax = pad->Charge();
438 distPrev = dist[indx];
439 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
440 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),
443 toBeRemoved.AddLast(pad);
444 fPreCluster->Pad(indx)->Release();
447 Int_t nRemove = toBeRemoved.GetEntriesFast();
448 for ( Int_t i = 0; i < nRemove; ++i )
450 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
453 } // if ( !cluster->IsSaturated() &&
457 AliDebug(2,"End of CheckPreClusterTwoCathodes=");
458 //StdoutToAliDebug(2,cluster->Print("full"));
463 //_____________________________________________________________________________
465 AliMUONClusterFinderPeakCOG::CheckOverlaps()
467 /// For debug only : check if some pixels overlap...
469 Int_t nPix = fPixArray->GetLast()+1;
472 for ( Int_t i = 0; i < nPix; ++i )
474 AliMUONPad* pixelI = Pixel(i);
475 AliMUONPad pi(dummy,dummy,dummy,dummy,
476 pixelI->Coord(0),pixelI->Coord(1),
477 pixelI->Size(0),pixelI->Size(1),0.0);
479 for ( Int_t j = i+1; j < nPix; ++j )
481 AliMUONPad* pixelJ = Pixel(j);
482 AliMUONPad pj(dummy,dummy,dummy,dummy,
483 pixelJ->Coord(0),pixelJ->Coord(1),
484 pixelJ->Size(0),pixelJ->Size(1),0.0);
487 if ( AliMUONPad::AreOverlapping(pi,pj,fgkDecreaseSize,area) )
489 AliInfo(Form("The following 2 pixels (%d and %d) overlap !",i,j));
491 StdoutToAliInfo(pixelI->Print();
492 cout << " Surface = " << pixelI->Size(0)*pixelI->Size(1)*4 << endl;
494 cout << " Surface = " << pixelJ->Size(0)*pixelJ->Size(1)*4 << endl;
495 cout << " Area surface = " << area.Dimensions().X()*area.Dimensions().Y()*4 << endl;
496 cout << "-------" << endl;
504 //_____________________________________________________________________________
505 void AliMUONClusterFinderPeakCOG::BuildPixArray(AliMUONCluster& cluster)
507 /// Build pixel array
509 Int_t npad = cluster.Multiplicity();
512 AliWarning("Got no pad at all ?!");
516 BuildPixArrayOneCathode(cluster);
518 // StdoutToAliDebug(2,cout << "End of BuildPixelArray:" << endl;
519 // fPixArray->Print(););
520 //CheckOverlaps();//FIXME : this is for debug only. Remove it.
523 //_____________________________________________________________________________
524 void AliMUONClusterFinderPeakCOG::BuildPixArrayOneCathode(AliMUONCluster& cluster)
526 /// Build the pixel array
528 // AliDebug(2,Form("cluster.Multiplicity=%d",cluster.Multiplicity()));
530 TVector2 dim = cluster.MinPadDimensions (-1, kFALSE);
531 Double_t width[2] = {dim.X(), dim.Y()}, xy0[2];
532 Int_t found[2] = {0}, mult = cluster.Multiplicity();
534 for ( Int_t i = 0; i < mult; ++i) {
535 AliMUONPad* pad = cluster.Pad(i);
536 for (Int_t j = 0; j < 2; ++j) {
537 if (found[j] == 0 && TMath::Abs(pad->Size(j)-width[j]) < fgkDistancePrecision) {
538 xy0[j] = pad->Coord(j);
542 if (found[0] && found[1]) break;
545 Double_t min[2], max[2];
546 Int_t cath0 = 0, cath1 = 1;
547 if (cluster.Multiplicity(0) == 0) cath0 = 1;
548 else if (cluster.Multiplicity(1) == 0) cath1 = 0;
550 TVector2 leftDown = cluster.Area(cath0).LeftDownCorner();
551 TVector2 rightUp = cluster.Area(cath0).RightUpCorner();
552 min[0] = leftDown.X();
553 min[1] = leftDown.Y();
554 max[0] = rightUp.X();
555 max[1] = rightUp.Y();
556 if (cath1 != cath0) {
557 leftDown = cluster.Area(cath1).LeftDownCorner();
558 rightUp = cluster.Area(cath1).RightUpCorner();
559 min[0] = TMath::Max (min[0], leftDown.X());
560 min[1] = TMath::Max (min[1], leftDown.Y());
561 max[0] = TMath::Min (max[0], rightUp.X());
562 max[1] = TMath::Min (max[1], rightUp.Y());
566 //width[0] /= 2; width[1] /= 2; // just for check
568 for (Int_t i = 0; i < 2; ++i) {
569 Double_t dist = (min[i] - xy0[i]) / width[i] / 2;
570 if (TMath::Abs(dist) < 1.e-6) dist = -1.e-6;
571 min[i] = xy0[i] + (TMath::Nint(dist-TMath::Sign(1.e-6,dist))
572 + TMath::Sign(0.5,dist)) * width[i] * 2;
573 nbins[i] = TMath::Nint ((max[i] - min[i]) / width[i] / 2);
574 if (nbins[i] == 0) ++nbins[i];
575 max[i] = min[i] + nbins[i] * width[i] * 2;
576 //cout << dist << " " << min[i] << " " << max[i] << " " << nbins[i] << endl;
580 TH2D *hist1 = new TH2D ("Grid", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
581 TH2D *hist2 = new TH2D ("Entries", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
582 TAxis *xaxis = hist1->GetXaxis();
583 TAxis *yaxis = hist1->GetYaxis();
586 for ( Int_t i = 0; i < mult; ++i) {
587 AliMUONPad* pad = cluster.Pad(i);
588 Int_t ix0 = xaxis->FindBin(pad->X());
589 Int_t iy0 = yaxis->FindBin(pad->Y());
590 PadOverHist(0, ix0, iy0, pad, hist1, hist2);
594 for (Int_t i = 1; i <= nbins[0]; ++i) {
595 Double_t x = xaxis->GetBinCenter(i);
596 for (Int_t j = 1; j <= nbins[1]; ++j) {
597 if (hist2->GetCellContent(i,j) < 0.1) continue;
598 if (cath0 != cath1) {
600 Double_t cont = hist2->GetCellContent(i,j);
601 if (cont < 999.) continue;
602 if (cont-Int_t(cont/1000.)*1000. < 0.5) continue;
604 //if (hist2->GetCellContent(i,j) < 1.1 && cluster.Multiplicity(0) &&
605 // cluster.Multiplicity(1)) continue;
606 Double_t y = yaxis->GetBinCenter(j);
607 Double_t charge = hist1->GetCellContent(i,j);
608 AliMUONPad* pixPtr = new AliMUONPad(x, y, width[0], width[1], charge);
609 fPixArray->Add(pixPtr);
613 if (fPixArray->GetEntriesFast() == 1) {
614 // Split pixel into 2
615 AliMUONPad* pixPtr = static_cast<AliMUONPad*> (fPixArray->UncheckedAt(0));
616 pixPtr->SetSize(0,width[0]/2.);
617 pixPtr->Shift(0,-width[0]/4.);
618 pixPtr = new AliMUONPad(pixPtr->X()+width[0], pixPtr->Y(), width[0]/2., width[1], pixPtr->Charge());
619 fPixArray->Add(pixPtr);
622 //fPixArray->Print();
627 //_____________________________________________________________________________
628 void AliMUONClusterFinderPeakCOG::PadOverHist(Int_t idir, Int_t ix0, Int_t iy0, AliMUONPad *pad,
629 TH2D *hist1, TH2D *hist2)
631 /// "Span" pad over histogram in the direction idir
633 TAxis *axis = idir == 0 ? hist1->GetXaxis() : hist1->GetYaxis();
634 Int_t nbins = axis->GetNbins(), cath = pad->Cathode();
635 Double_t bin = axis->GetBinWidth(1), amask = TMath::Power(1000.,cath*1.);
637 Int_t nbinPad = (Int_t)(pad->Size(idir)/bin*2+fgkDistancePrecision) + 1; // number of bins covered by pad
639 for (Int_t i = 0; i < nbinPad; ++i) {
640 Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
641 if (ixy > nbins) break;
642 Double_t lowEdge = axis->GetBinLowEdge(ixy);
643 if (lowEdge + fgkDistancePrecision > pad->Coord(idir) + pad->Size(idir)) break;
644 if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
647 Double_t cont = pad->Charge();
648 if (hist2->GetCellContent(ix0, ixy) > 0.1)
649 cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont)
650 + TMath::Min (TMath::Max(hist1->GetCellContent(ix0, ixy),cont)*0.1, 10.);
651 hist1->SetCellContent(ix0, ixy, cont);
652 hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
656 for (Int_t i = -1; i > -nbinPad; --i) {
657 Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
659 Double_t upEdge = axis->GetBinUpEdge(ixy);
660 if (upEdge - fgkDistancePrecision < pad->Coord(idir) - pad->Size(idir)) break;
661 if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
664 Double_t cont = pad->Charge();
665 if (hist2->GetCellContent(ix0, ixy) > 0.1)
666 cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont)
667 + TMath::Min (TMath::Max(hist1->GetCellContent(ix0, ixy),cont)*0.1,10.);
668 hist1->SetCellContent(ix0, ixy, cont);
669 hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
674 //_____________________________________________________________________________
675 Int_t AliMUONClusterFinderPeakCOG::FindLocalMaxima(TObjArray *pixArray, Int_t *localMax, Double_t *maxVal)
677 /// Find local maxima in pixel space
679 AliDebug(1,Form("nPix=%d",pixArray->GetLast()+1));
682 //delete ((TH2D*) gROOT->FindObject("anode"));
683 //if (pixArray == fPixArray) hist = (TH2D*) gROOT->FindObject("anode");
684 //else { hist = (TH2D*) gROOT->FindObject("anode1"); cout << hist << endl; }
685 //if (hist) hist->Delete();
688 Double_t xylim[4] = {999, 999, 999, 999};
690 Int_t nPix = pixArray->GetEntriesFast();
691 AliMUONPad *pixPtr = 0;
692 for (Int_t ipix = 0; ipix < nPix; ++ipix) {
693 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
694 for (Int_t i = 0; i < 4; ++i)
695 xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
697 for (Int_t i = 0; i < 4; ++i) xylim[i] -= pixPtr->Size(i/2);
699 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
700 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
701 if (pixArray == fPixArray) fHistAnode = new TH2D("anode","anode",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
702 else fHistAnode = new TH2D("anode1","anode1",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
703 for (Int_t ipix = 0; ipix < nPix; ++ipix) {
704 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
705 fHistAnode->Fill(pixPtr->Coord(0), pixPtr->Coord(1), pixPtr->Charge());
707 // if (fDraw && pixArray == fPixArray) fDraw->DrawHist("c2", hist);
709 Int_t nMax = 0, indx, nxy = ny * nx;
710 Int_t *isLocalMax = new Int_t[nxy];
711 for (Int_t i = 0; i < nxy; ++i) isLocalMax[i] = 0;
713 for (Int_t i = 1; i <= ny; ++i) {
715 for (Int_t j = 1; j <= nx; ++j) {
716 if (fHistAnode->GetCellContent(j,i) < 0.5) continue;
717 //if (isLocalMax[indx+j-1] < 0) continue;
718 if (isLocalMax[indx+j-1] != 0) continue;
719 FlagLocalMax(fHistAnode, i, j, isLocalMax);
723 for (Int_t i = 1; i <= ny; ++i) {
725 for (Int_t j = 1; j <= nx; ++j) {
726 if (isLocalMax[indx+j-1] > 0) {
727 localMax[nMax] = indx + j - 1;
728 maxVal[nMax++] = fHistAnode->GetCellContent(j,i);
729 if (nMax > 99) break;
733 AliError(" Too many local maxima !!!");
737 if (fDebug) cout << " Local max: " << nMax << endl;
738 delete [] isLocalMax;
742 //_____________________________________________________________________________
743 void AliMUONClusterFinderPeakCOG::FlagLocalMax(TH2D *hist, Int_t i, Int_t j, Int_t *isLocalMax)
745 /// Flag pixels (whether or not local maxima)
747 Int_t nx = hist->GetNbinsX();
748 Int_t ny = hist->GetNbinsY();
749 Int_t cont = TMath::Nint (hist->GetCellContent(j,i));
750 Int_t cont1 = 0, indx = (i-1)*nx+j-1, indx1 = 0, indx2 = 0;
752 Int_t ie = i + 2, je = j + 2;
753 for (Int_t i1 = i-1; i1 < ie; ++i1) {
754 if (i1 < 1 || i1 > ny) continue;
755 indx1 = (i1 - 1) * nx;
756 for (Int_t j1 = j-1; j1 < je; ++j1) {
757 if (j1 < 1 || j1 > nx) continue;
758 if (i == i1 && j == j1) continue;
759 indx2 = indx1 + j1 - 1;
760 cont1 = TMath::Nint (hist->GetCellContent(j1,i1));
761 if (cont < cont1) { isLocalMax[indx] = -1; return; }
762 else if (cont > cont1) isLocalMax[indx2] = -1;
763 else { // the same charge
764 isLocalMax[indx] = 1;
765 if (isLocalMax[indx2] == 0) {
766 FlagLocalMax(hist, i1, j1, isLocalMax);
767 if (isLocalMax[indx2] < 0) { isLocalMax[indx] = -1; return; }
768 else isLocalMax[indx2] = -1;
773 isLocalMax[indx] = 1; // local maximum
776 //_____________________________________________________________________________
777 void AliMUONClusterFinderPeakCOG::FindCluster(AliMUONCluster& cluster,
778 Int_t *localMax, Int_t iMax)
780 /// Find pixel cluster around local maximum \a iMax and pick up pads
781 /// overlapping with it
783 //TH2D *hist = (TH2D*) gROOT->FindObject("anode");
785 TCanvas* c = new TCanvas("Anode","Anode",800,600);
787 hist->Draw("lego1Fb"); // debug
792 Int_t nx = fHistAnode->GetNbinsX();
793 //Int_t ny = hist->GetNbinsY();
794 Int_t ic = localMax[iMax] / nx + 1;
795 Int_t jc = localMax[iMax] % nx + 1;
797 // Get min pad dimensions for the precluster
799 if (cluster.Multiplicity(0) == 0 || cluster.Multiplicity(1) == 0) nSides = 1;
800 TVector2 dim0 = cluster.MinPadDimensions(0, -1, kFALSE);
801 TVector2 dim1 = cluster.MinPadDimensions(1, -1, kFALSE);
802 //Double_t width[2][2] = {{dim0.X(), dim0.Y()},{dim1.X(),dim1.Y()}};
803 Int_t nonb[2] = {1, 0}; // coordinate index vs cathode
804 if (nSides == 1 || dim0.X() < dim1.X() - fgkDistancePrecision) {
808 Bool_t samex = kFALSE, samey = kFALSE;
809 if (TMath::Abs(dim0.X()-dim1.X()) < fgkDistancePrecision) samex = kTRUE; // the same X pad size on both planes
810 if (TMath::Abs(dim0.Y()-dim1.Y()) < fgkDistancePrecision) samey = kTRUE; // the same Y pad size on both planes
812 // Drop all pixels from the array - pick up only the ones from the cluster
813 //fPixArray->Delete();
815 Double_t wx = fHistAnode->GetXaxis()->GetBinWidth(1)/2;
816 Double_t wy = fHistAnode->GetYaxis()->GetBinWidth(1)/2;
817 Double_t yc = fHistAnode->GetYaxis()->GetBinCenter(ic);
818 Double_t xc = fHistAnode->GetXaxis()->GetBinCenter(jc);
819 Double_t cont = fHistAnode->GetCellContent(jc,ic);
820 AliMUONPad* pixPtr = new AliMUONPad (xc, yc, wx, wy, cont);
821 if (fDebug) pixPtr->Print("full");
823 Int_t npad = cluster.Multiplicity();
825 // Pick up pads which overlap with the maximum pixel and find pads with the max signal
826 Double_t qMax[2] = {0};
827 AliMUONPad *matrix[2][3] = {{0x0,0x0,0x0},{0x0,0x0,0x0}};
828 for (Int_t j = 0; j < npad; ++j)
830 AliMUONPad* pad = cluster.Pad(j);
831 if ( Overlap(*pad,*pixPtr) )
833 if (fDebug) { cout << j << " "; pad->Print("full"); }
834 if (pad->Charge() > qMax[pad->Cathode()]) {
835 qMax[pad->Cathode()] = pad->Charge();
836 matrix[pad->Cathode()][1] = pad;
837 if (nSides == 1) matrix[!pad->Cathode()][1] = pad;
841 //if (nSides == 2 && (matrix[0][1] == 0x0 || matrix[1][1] == 0x0)) return; // ???
843 // Find neighbours of maxima to have 3 pads per direction (if possible)
844 for (Int_t j = 0; j < npad; ++j)
846 AliMUONPad* pad = cluster.Pad(j);
847 Int_t cath = pad->Cathode();
848 if (pad == matrix[cath][1]) continue;
849 Int_t nLoops = 3 - nSides;
851 for (Int_t k = 0; k < nLoops; ++k) {
853 if (k) cath1 = !cath;
855 // Check the coordinate corresponding to the cathode (bending or non-bending case)
856 Double_t dist = pad->Coord(nonb[cath1]) - matrix[cath][1]->Coord(nonb[cath1]);
857 Double_t dir = TMath::Sign (1., dist);
858 dist = TMath::Abs(dist) - pad->Size(nonb[cath1]) - matrix[cath][1]->Size(nonb[cath1]);
860 if (TMath::Abs(dist) < fgkDistancePrecision) {
861 // Check the other coordinate
862 dist = pad->Coord(!nonb[cath1]) - matrix[cath1][1]->Coord(!nonb[cath1]);
863 if (TMath::Abs(dist) >
864 TMath::Max(pad->Size(!nonb[cath1]), matrix[cath1][1]->Size(!nonb[cath1])) - fgkDistancePrecision) break;
865 Int_t idir = TMath::Nint (dir);
866 if (matrix[cath1][1+idir] == 0x0) matrix[cath1][1+idir] = pad;
867 else if (pad->Charge() > matrix[cath1][1+idir]->Charge()) matrix[cath1][1+idir] = pad; // diff. segmentation
868 //cout << pad->Coord(nonb[cath1]) << " " << pad->Coord(!nonb[cath1]) << " " << pad->Size(nonb[cath1]) << " " << pad->Size(!nonb[cath1]) << " " << pad->Charge() << endl ;
874 Double_t coord[2] = {0.}, qAver = 0.;
875 for (Int_t i = 0; i < 2; ++i) {
877 Double_t coordQ = 0.;
878 Int_t cath = matrix[i][1]->Cathode();
879 if (i && nSides == 1) cath = !cath;
880 for (Int_t j = 0; j < 3; ++j) {
881 if (matrix[i][j] == 0x0) continue;
882 Double_t dq = matrix[i][j]->Charge();
884 coordQ += dq * matrix[i][j]->Coord(nonb[cath]);
885 //coordQ += (matrix[i][j]->Charge() * matrix[i][j]->Coord(nonb[cath]));
887 coord[cath] = coordQ / q;
888 qAver = TMath::Max (qAver, q);
891 //qAver = TMath::Sqrt(qAver);
895 AliMUONCluster* cluster1 = new AliMUONCluster(cluster);
897 cluster1->SetCharge(qAver,qAver);
899 cluster1->SetPosition(TVector2(coord[1],coord[0]),TVector2(0.,0.));
901 cluster1->SetPosition(TVector2(coord[0],coord[1]),TVector2(0.,0.));
903 cluster1->SetChi2(0.);
905 // FIXME: we miss some information in this cluster, as compared to
906 // the original AddRawCluster code.
908 AliDebug(2,Form("Adding RawCluster detElemId %4d mult %2d charge %5d (xl,yl)=(%9.6g,%9.6g)",
909 fDetElemId,cluster1->Multiplicity(),(Int_t)cluster1->Charge(),
910 cluster1->Position().X(),cluster1->Position().Y()));
912 fClusterList.Add(cluster1);
916 //_____________________________________________________________________________
917 AliMUONClusterFinderPeakCOG&
918 AliMUONClusterFinderPeakCOG::operator=(const AliMUONClusterFinderPeakCOG& rhs)
920 /// Protected assignement operator
922 if (this == &rhs) return *this;
924 AliFatal("Not implemented.");
929 //_____________________________________________________________________________
930 void AliMUONClusterFinderPeakCOG::PadsInXandY(AliMUONCluster& cluster,
931 Int_t &nInX, Int_t &nInY) const
933 /// Find number of pads in X and Y-directions (excluding virtual ones and
936 //Int_t statusToTest = 1;
937 Int_t statusToTest = fgkUseForFit;
939 //if ( nInX < 0 ) statusToTest = 0;
940 if ( nInX < 0 ) statusToTest = fgkZero;
942 Bool_t mustMatch(kTRUE);
944 AliMpIntPair cn = cluster.NofPads(statusToTest,mustMatch);
946 nInX = cn.GetFirst();
947 nInY = cn.GetSecond();
950 //_____________________________________________________________________________
951 void AliMUONClusterFinderPeakCOG::RemovePixel(Int_t i)
953 /// Remove pixel at index i
954 AliMUONPad* pixPtr = Pixel(i);
955 fPixArray->RemoveAt(i);
959 //_____________________________________________________________________________
961 AliMUONClusterFinderPeakCOG::Pixel(Int_t i) const
963 /// Returns pixel at index i
964 return static_cast<AliMUONPad*>(fPixArray->UncheckedAt(i));
967 //_____________________________________________________________________________
969 AliMUONClusterFinderPeakCOG::Print(Option_t* what) const
974 if ( swhat.Contains("precluster") )
976 if ( fPreCluster) fPreCluster->Print();