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 AliMUONClusterFinderMLEM
21 /// Clusterizer class based on the Expectation-Maximization algorithm
23 /// Pre-clustering is handled by AliMUONPreClusterFinder
24 /// From a precluster a pixel array is built, and from this pixel array
25 /// a list of clusters is output (using AliMUONClusterSplitterMLEM).
27 /// \author Laurent Aphecetche (for the "new" C++ structure) and
28 /// Alexander Zinchenko, JINR Dubna, for the hardcore of it ;-)
29 //-----------------------------------------------------------------------------
31 #include "AliMUONClusterFinderMLEM.h"
33 #include "AliMUONCluster.h"
34 #include "AliMUONClusterSplitterMLEM.h"
35 #include "AliMUONVDigit.h"
36 #include "AliMUONPad.h"
37 #include "AliMUONPreClusterFinder.h"
39 #include "AliMpVPadIterator.h"
40 #include "AliMpVSegmentation.h"
41 #include "AliRunLoader.h"
42 #include "AliMUONVDigitStore.h"
43 #include <Riostream.h>
48 #include "AliCodeTimer.h"
51 ClassImp(AliMUONClusterFinderMLEM)
54 const Double_t AliMUONClusterFinderMLEM::fgkZeroSuppression = 6; // average zero suppression value
55 //const Double_t AliMUONClusterFinderMLEM::fgkDistancePrecision = 1e-6; // (cm) used to check overlaps and so on
56 const Double_t AliMUONClusterFinderMLEM::fgkDistancePrecision = 1e-3; // (cm) used to check overlaps and so on
57 const TVector2 AliMUONClusterFinderMLEM::fgkIncreaseSize(-AliMUONClusterFinderMLEM::fgkDistancePrecision,-AliMUONClusterFinderMLEM::fgkDistancePrecision);
58 const TVector2 AliMUONClusterFinderMLEM::fgkDecreaseSize(AliMUONClusterFinderMLEM::fgkDistancePrecision,AliMUONClusterFinderMLEM::fgkDistancePrecision);
60 // Status flags for pads
61 const Int_t AliMUONClusterFinderMLEM::fgkZero = 0x0; ///< pad "basic" state
62 const Int_t AliMUONClusterFinderMLEM::fgkMustKeep = 0x1; ///< do not kill (for pixels)
63 const Int_t AliMUONClusterFinderMLEM::fgkUseForFit = 0x10; ///< should be used for fit
64 const Int_t AliMUONClusterFinderMLEM::fgkOver = 0x100; ///< processing is over
65 const Int_t AliMUONClusterFinderMLEM::fgkModified = 0x1000; ///< modified pad charge
66 const Int_t AliMUONClusterFinderMLEM::fgkCoupled = 0x10000; ///< coupled pad
68 //_____________________________________________________________________________
69 AliMUONClusterFinderMLEM::AliMUONClusterFinderMLEM(Bool_t plot, AliMUONVClusterFinder* clusterFinder)
70 : AliMUONVClusterFinder(),
71 fPreClusterFinder(clusterFinder),
79 fPixArray(new TObjArray(20)),
88 fkSegmentation[1] = fkSegmentation[0] = 0x0;
90 if (fPlot) fDebug = 1;
93 //_____________________________________________________________________________
94 AliMUONClusterFinderMLEM::~AliMUONClusterFinderMLEM()
97 delete fPixArray; fPixArray = 0;
99 delete fPreClusterFinder;
101 AliInfo(Form("Total clusters %d AddVirtualPad needed %d",
102 fNClusters,fNAddVirtualPads));
105 //_____________________________________________________________________________
107 AliMUONClusterFinderMLEM::Prepare(Int_t detElemId,
108 TClonesArray* pads[2],
109 const AliMpArea& area,
110 const AliMpVSegmentation* seg[2])
112 /// Prepare for clustering
113 // AliCodeTimerAuto("")
115 for ( Int_t i = 0; i < 2; ++i )
117 fkSegmentation[i] = seg[i];
120 // Find out the DetElemId
121 fDetElemId = detElemId;
124 fSplitter = new AliMUONClusterSplitterMLEM(fDetElemId,fPixArray);
125 fSplitter->SetDebug(fDebug);
127 // find out current event number, and reset the cluster number
128 AliRunLoader *runLoader = AliRunLoader::Instance();
129 fEventNumber = runLoader ? runLoader->GetEventNumber() : 0;
131 fClusterList.Delete();
133 AliDebug(3,Form("EVT %d DE %d",fEventNumber,fDetElemId));
135 if ( fPreClusterFinder->NeedSegmentation() )
137 return fPreClusterFinder->Prepare(detElemId,pads,area,seg);
141 return fPreClusterFinder->Prepare(detElemId,pads,area);
145 //_____________________________________________________________________________
147 AliMUONClusterFinderMLEM::NextCluster()
149 /// Return next cluster
150 // AliCodeTimerAuto("")
152 // if the list of clusters is not void, pick one from there
153 TObject* o = fClusterList.At(++fClusterNumber);
154 if ( o != 0x0 ) return static_cast<AliMUONCluster*>(o);
156 //FIXME : at this point, must check whether we've used all the digits
157 //from precluster : if not, let the preclustering know about those unused
158 //digits, so it can reuse them
160 // if the cluster list is exhausted, we need to go to the next
161 // pre-cluster and treat it
163 fPreCluster = fPreClusterFinder->NextCluster();
165 fClusterList.Delete(); // reset the list of clusters for this pre-cluster
166 fClusterNumber = -1; //AZ
176 // WorkOnPreCluster may have used only part of the pads, so we check that
177 // now, and let the unused pads be reused by the preclustering...
179 Int_t mult = fPreCluster->Multiplicity();
180 for ( Int_t i = 0; i < mult; ++i )
182 AliMUONPad* pad = fPreCluster->Pad(i);
183 if ( !pad->IsUsed() )
185 fPreClusterFinder->UsePad(*pad);
189 return NextCluster();
192 //_____________________________________________________________________________
194 AliMUONClusterFinderMLEM::WorkOnPreCluster()
196 /// Starting from a precluster, builds a pixel array, and then
197 /// extract clusters from this array
199 // AliCodeTimerAuto("")
202 cout << " *** Event # " << fEventNumber
203 << " det. elem.: " << fDetElemId << endl;
204 for (Int_t j = 0; j < fPreCluster->Multiplicity(); ++j) {
205 AliMUONPad* pad = fPreCluster->Pad(j);
206 printf(" bbb %3d %1d %8.4f %8.4f %8.4f %8.4f %6.1f %3d %3d %2d %1d %1d \n",
207 j, pad->Cathode(), pad->Coord(0), pad->Coord(1), pad->DX()*2, pad->DY()*2,
208 pad->Charge(), pad->Ix(), pad->Iy(), pad->Status(), pad->IsReal(), pad->IsSaturated());
212 AliMUONCluster* cluster = CheckPrecluster(*fPreCluster);
213 if (!cluster) return kFALSE;
215 BuildPixArray(*cluster);
217 if ( fPixArray->GetLast() < 0 )
219 AliDebug(1,"No pixel for the above cluster");
224 // Use MLEM for cluster finder
225 Int_t nMax = 1, localMax[100], maxPos[100];
226 Double_t maxVal[100];
228 Int_t iSimple = 0, nInX = -1, nInY;
230 PadsInXandY(*cluster,nInX, nInY);
232 if (nInX < 4 && nInY < 4)
234 iSimple = 1; // simple cluster
238 nMax = FindLocalMaxima(fPixArray, localMax, maxVal); // for small clusters just to tag pixels
240 if (cluster->Multiplicity() <= 50) nMax = 1; // for small clusters
241 if (nMax > 1) TMath::Sort(nMax, maxVal, maxPos, kTRUE); // in descending order
245 for (Int_t i = 0; i < nMax; ++i)
249 FindCluster(*cluster,localMax, maxPos[i]);
252 MainLoop(*cluster,iSimple);
256 Int_t mult = cluster->Multiplicity();
257 for (Int_t j = 0; j < mult; ++j)
259 AliMUONPad* pad = cluster->Pad(j);
260 //if ( pad->Status() == 0 ) continue; // pad charge was not modified
261 if ( pad->Status() != fgkOver ) continue; // pad was not used
263 pad->SetStatus(fgkZero);
264 pad->RevertCharge(); // use backup charge value
267 } // for (Int_t i=0; i<nMax;
270 fHistMlem = fHistAnode = 0x0;
275 //_____________________________________________________________________________
277 AliMUONClusterFinderMLEM::Overlap(const AliMUONPad& pad, const AliMUONPad& pixel)
279 /// Check if the pad and the pixel overlaps
281 // make a fake pad from the pixel
282 AliMUONPad tmp(pad.DetElemId(),pad.Cathode(),pad.Ix(),pad.Iy(),
283 pixel.Coord(0),pixel.Coord(1),
284 pixel.Size(0),pixel.Size(1),0);
286 return AliMUONPad::AreOverlapping(pad,tmp,fgkDecreaseSize);
289 //_____________________________________________________________________________
291 AliMUONClusterFinderMLEM::CheckPrecluster(const AliMUONCluster& origCluster)
293 /// Check precluster in order to attempt to simplify it (mostly for
294 /// two-cathode preclusters)
298 // Disregard small clusters (leftovers from splitting or noise)
299 if ((origCluster.Multiplicity()==1 || origCluster.Multiplicity()==2) &&
300 origCluster.Charge(0)+origCluster.Charge(1) < 10)
305 AliMUONCluster* cluster = new AliMUONCluster(origCluster);
307 AliDebug(2,"Start of CheckPreCluster=");
308 //StdoutToAliDebug(2,cluster->Print("full"));
310 AliMUONCluster* rv(0x0);
312 if (cluster->Multiplicity(0) && cluster->Multiplicity(1))
314 rv = CheckPreclusterTwoCathodes(cluster);
323 //_____________________________________________________________________________
325 AliMUONClusterFinderMLEM::CheckPreclusterTwoCathodes(AliMUONCluster* cluster)
327 /// Check two-cathode cluster
329 Int_t npad = cluster->Multiplicity();
330 Int_t* flags = new Int_t[npad];
331 for (Int_t j = 0; j < npad; ++j) flags[j] = 0;
333 // Check pad overlaps
334 for ( Int_t i = 0; i < npad; ++i)
336 AliMUONPad* padi = cluster->Pad(i);
337 if ( padi->Cathode() != 0 ) continue;
338 for (Int_t j = i+1; j < npad; ++j)
340 AliMUONPad* padj = cluster->Pad(j);
341 if ( padj->Cathode() != 1 ) continue;
342 if ( !AliMUONPad::AreOverlapping(*padi,*padj,fgkDecreaseSize) ) continue;
343 flags[i] = flags[j] = 1; // mark overlapped pads
347 // Check if all pads overlap
349 for (Int_t i = 0; i < npad; ++i)
351 if (!flags[i]) ++nFlags;
356 // not all pads overlap.
357 if (fDebug) cout << " nFlags: " << nFlags << endl;
358 TObjArray toBeRemoved;
359 for (Int_t i = 0; i < npad; ++i)
361 AliMUONPad* pad = cluster->Pad(i);
362 if (flags[i]) continue;
363 Int_t cath = pad->Cathode();
364 Int_t cath1 = TMath::Even(cath);
365 // Check for edge effect (missing pads on the _other_ cathode)
367 fkSegmentation[cath1]->PadByPosition(pad->Position().X(),
368 pad->Position().Y(),kFALSE);
369 if (!mpPad.IsValid()) continue;
370 //if (nFlags == 1 && pad->Charge() < fgkZeroSuppression * 3) continue;
371 if (nFlags == 1 && pad->Charge() < 20) continue;
372 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
373 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),pad->Charge()));
374 toBeRemoved.AddLast(pad);
375 fPreCluster->Pad(i)->Release();
377 Int_t nRemove = toBeRemoved.GetEntriesFast();
378 for ( Int_t i = 0; i < nRemove; ++i )
380 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
384 // Check correlations of cathode charges
385 if ( !cluster->IsSaturated() && cluster->ChargeAsymmetry() > 1 )
388 Int_t cathode = cluster->MaxRawChargeCathode();
394 // get min and max pad charges on the cathode opposite to the
395 // max pad (given by MaxRawChargeCathode())
397 Int_t mult = cluster->Multiplicity();
398 for ( Int_t i = 0; i < mult; ++i )
400 AliMUONPad* pad = cluster->Pad(i);
401 if ( pad->Cathode() != cathode || !pad->IsReal() )
403 // only consider pads in the opposite cathode, and
404 // only consider real pads (i.e. exclude the virtual ones)
407 if ( pad->Charge() < cmin )
409 cmin = pad->Charge();
416 else if ( pad->Charge() > cmax )
418 cmax = pad->Charge();
422 AliDebug(2,Form("Pad imin,imax %d,%d cmin,cmax %e,%e",
423 imin,imax,cmin,cmax));
425 // arrange pads according to their distance to the max, normalized
427 Double_t* dist = new Double_t[mult];
432 AliMUONPad* padmax = cluster->Pad(imax);
434 for ( Int_t i = 0; i < mult; ++i )
437 if ( i == imax) continue;
438 AliMUONPad* pad = cluster->Pad(i);
439 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
440 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
441 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
442 dist[i] = TMath::Sqrt(dx*dx+dy*dy);
445 dmin = dist[i] + 1E-3; // distance to the pad with minimum charge
451 TMath::Sort(mult,dist,flags,kFALSE); // in ascending order
452 Double_t xmax(-1), distPrev(999);
453 TObjArray toBeRemoved;
455 for ( Int_t i = 0; i < mult; ++i )
457 Int_t indx = flags[i];
458 AliMUONPad* pad = cluster->Pad(indx);
459 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
460 if ( dist[indx] > dmin )
462 // farther than the minimum pad
463 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
464 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
467 if (dx >= 0 && dy >= 0) continue;
468 if (TMath::Abs(dx) > TMath::Abs(dy) && dx >= 0) continue;
469 if (TMath::Abs(dy) > TMath::Abs(dx) && dy >= 0) continue;
471 if (dist[indx] > distPrev + 1) break; // overstepping empty pads
472 if ( pad->Charge() <= cmax || TMath::Abs(dist[indx]-xmax) < 1E-3 )
475 if (TMath::Abs(dist[indx]-xmax) < 1.e-3)
477 cmax = TMath::Max(pad->Charge(),cmax);
481 cmax = pad->Charge();
484 distPrev = dist[indx];
485 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
486 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),
489 toBeRemoved.AddLast(pad);
490 fPreCluster->Pad(indx)->Release();
493 Int_t nRemove = toBeRemoved.GetEntriesFast();
494 for ( Int_t i = 0; i < nRemove; ++i )
496 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
503 AliDebug(2,"End of CheckPreClusterTwoCathodes=");
504 //StdoutToAliDebug(2,cluster->Print("full"));
509 //_____________________________________________________________________________
511 AliMUONClusterFinderMLEM::CheckOverlaps()
513 /// For debug only : check if some pixels overlap...
515 Int_t nPix = fPixArray->GetLast()+1;
518 for ( Int_t i = 0; i < nPix; ++i )
520 AliMUONPad* pixelI = Pixel(i);
521 AliMUONPad pi(dummy,dummy,dummy,dummy,
522 pixelI->Coord(0),pixelI->Coord(1),
523 pixelI->Size(0),pixelI->Size(1),0.0);
525 for ( Int_t j = i+1; j < nPix; ++j )
527 AliMUONPad* pixelJ = Pixel(j);
528 AliMUONPad pj(dummy,dummy,dummy,dummy,
529 pixelJ->Coord(0),pixelJ->Coord(1),
530 pixelJ->Size(0),pixelJ->Size(1),0.0);
533 if ( AliMUONPad::AreOverlapping(pi,pj,fgkDecreaseSize,area) )
535 AliInfo(Form("The following 2 pixels (%d and %d) overlap !",i,j));
537 StdoutToAliInfo(pixelI->Print();
538 cout << " Surface = " << pixelI->Size(0)*pixelI->Size(1)*4 << endl;
540 cout << " Surface = " << pixelJ->Size(0)*pixelJ->Size(1)*4 << endl;
541 cout << " Area surface = " << area.GetDimensionX()*area.
\14\14*4 << endl;
542 cout << "-------" << endl;
550 //_____________________________________________________________________________
551 void AliMUONClusterFinderMLEM::BuildPixArray(AliMUONCluster& cluster)
553 /// Build pixel array for MLEM method
555 Int_t npad = cluster.Multiplicity();
558 AliWarning("Got no pad at all ?!");
562 BuildPixArrayOneCathode(cluster);
564 Int_t nPix = fPixArray->GetLast()+1;
566 // AliDebug(2,Form("nPix after BuildPixArray=%d",nPix));
570 // AliDebug(2,Form("Will trim number of pixels to number of pads"));
572 // Too many pixels - sort and remove pixels with the lowest signal
574 for ( Int_t i = npad; i < nPix; ++i )
578 fPixArray->Compress();
579 } // if (nPix > npad)
581 // StdoutToAliDebug(2,cout << "End of BuildPixelArray:" << endl;
582 // fPixArray->Print(););
583 //CheckOverlaps();//FIXME : this is for debug only. Remove it.
586 //_____________________________________________________________________________
587 void AliMUONClusterFinderMLEM::BuildPixArrayOneCathode(AliMUONCluster& cluster)
589 /// Build the pixel array
591 // AliDebug(2,Form("cluster.Multiplicity=%d",cluster.Multiplicity()));
593 TVector2 dim = cluster.MinPadDimensions (-1, kFALSE);
594 Double_t width[2] = {dim.X(), dim.Y()}, xy0[2]={99999,99999};
595 Int_t found[2] = {0,0}, mult = cluster.Multiplicity();
597 for ( Int_t i = 0; i < mult; ++i) {
598 AliMUONPad* pad = cluster.Pad(i);
599 for (Int_t j = 0; j < 2; ++j) {
600 if (found[j] == 0 && TMath::Abs(pad->Size(j)-width[j]) < fgkDistancePrecision) {
601 xy0[j] = pad->Coord(j);
605 if (found[0] && found[1]) break;
608 Double_t min[2], max[2];
609 Int_t cath0 = 0, cath1 = 1;
610 if (cluster.Multiplicity(0) == 0) cath0 = 1;
611 else if (cluster.Multiplicity(1) == 0) cath1 = 0;
614 Double_t leftDownX, leftDownY;
615 cluster.Area(cath0).LeftDownCorner(leftDownX, leftDownY);
616 Double_t rightUpX, rightUpY;
617 cluster.Area(cath0).RightUpCorner(rightUpX, rightUpY);
622 if (cath1 != cath0) {
623 cluster.Area(cath1).LeftDownCorner(leftDownX, leftDownY);
624 cluster.Area(cath1).RightUpCorner(rightUpX, rightUpY);
625 min[0] = TMath::Max (min[0], leftDownX);
626 min[1] = TMath::Max (min[1], leftDownY);
627 max[0] = TMath::Min (max[0], rightUpX);
628 max[1] = TMath::Min (max[1], rightUpY);
632 //width[0] /= 2; width[1] /= 2; // just for check
633 Int_t nbins[2]={0,0};
634 for (Int_t i = 0; i < 2; ++i) {
635 Double_t dist = (min[i] - xy0[i]) / width[i] / 2;
636 if (TMath::Abs(dist) < 1.e-6) dist = -1.e-6;
637 min[i] = xy0[i] + (TMath::Nint(dist-TMath::Sign(1.e-6,dist))
638 + TMath::Sign(0.5,dist)) * width[i] * 2;
639 nbins[i] = TMath::Nint ((max[i] - min[i]) / width[i] / 2);
640 if (nbins[i] == 0) ++nbins[i];
641 max[i] = min[i] + nbins[i] * width[i] * 2;
642 //cout << dist << " " << min[i] << " " << max[i] << " " << nbins[i] << endl;
646 TH2D *hist1 = new TH2D ("Grid", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
647 TH2D *hist2 = new TH2D ("Entries", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
648 TAxis *xaxis = hist1->GetXaxis();
649 TAxis *yaxis = hist1->GetYaxis();
652 for ( Int_t i = 0; i < mult; ++i) {
653 AliMUONPad* pad = cluster.Pad(i);
654 Int_t ix0 = xaxis->FindBin(pad->X());
655 Int_t iy0 = yaxis->FindBin(pad->Y());
656 PadOverHist(0, ix0, iy0, pad, hist1, hist2);
660 for (Int_t i = 1; i <= nbins[0]; ++i) {
661 Double_t x = xaxis->GetBinCenter(i);
662 for (Int_t j = 1; j <= nbins[1]; ++j) {
663 if (hist2->GetCellContent(i,j) < 0.1) continue;
664 //if (hist2->GetCellContent(i,j) < 1.1 && cluster.Multiplicity(0) &&
665 // cluster.Multiplicity(1)) continue;
666 if (cath0 != cath1) {
668 Double_t cont = hist2->GetCellContent(i,j);
669 if (cont < 999.) continue;
670 if (cont-Int_t(cont/1000.)*1000. < 0.5) continue;
672 Double_t y = yaxis->GetBinCenter(j);
673 Double_t charge = hist1->GetCellContent(i,j);
674 AliMUONPad* pixPtr = new AliMUONPad(x, y, width[0], width[1], charge);
675 fPixArray->Add(pixPtr);
679 if (fPixArray->GetEntriesFast() == 1) {
680 // Split pixel into 2
681 AliMUONPad* pixPtr = static_cast<AliMUONPad*> (fPixArray->UncheckedAt(0));
682 pixPtr->SetSize(0,width[0]/2.);
683 pixPtr->Shift(0,-width[0]/4.);
684 pixPtr = new AliMUONPad(pixPtr->X()+width[0], pixPtr->Y(), width[0]/2., width[1], pixPtr->Charge());
685 fPixArray->Add(pixPtr);
688 //fPixArray->Print();
693 //_____________________________________________________________________________
694 void AliMUONClusterFinderMLEM::PadOverHist(Int_t idir, Int_t ix0, Int_t iy0, AliMUONPad *pad,
695 TH2D *hist1, TH2D *hist2)
697 /// "Span" pad over histogram in the direction idir
699 TAxis *axis = idir == 0 ? hist1->GetXaxis() : hist1->GetYaxis();
700 Int_t nbins = axis->GetNbins(), cath = pad->Cathode();
701 Double_t bin = axis->GetBinWidth(1), amask = TMath::Power(1000.,cath*1.);
703 Int_t nbinPad = (Int_t)(pad->Size(idir)/bin*2+fgkDistancePrecision) + 1; // number of bins covered by pad
705 for (Int_t i = 0; i < nbinPad; ++i) {
706 Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
707 if (ixy > nbins) break;
708 Double_t lowEdge = axis->GetBinLowEdge(ixy);
709 if (lowEdge + fgkDistancePrecision > pad->Coord(idir) + pad->Size(idir)) break;
710 if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
713 Double_t cont = pad->Charge();
714 if (hist2->GetCellContent(ix0, ixy) > 0.1)
715 cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont);
716 hist1->SetCellContent(ix0, ixy, cont);
717 //hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+1);
718 hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
722 for (Int_t i = -1; i > -nbinPad; --i) {
723 Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
725 Double_t upEdge = axis->GetBinUpEdge(ixy);
726 if (upEdge - fgkDistancePrecision < pad->Coord(idir) - pad->Size(idir)) break;
727 if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
730 Double_t cont = pad->Charge();
731 if (hist2->GetCellContent(ix0, ixy) > 0.1)
732 cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont);
733 hist1->SetCellContent(ix0, ixy, cont);
734 //hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+1);
735 hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
740 //_____________________________________________________________________________
742 AliMUONClusterFinderMLEM::Plot(const char* basename)
744 /// Make a plot and save it as png
749 TCanvas* c = new TCanvas("MLEM","MLEM",800,600);
754 c->Print(Form("%s.EVT%d.DE%d.CLU%d.png",basename,fEventNumber,
755 fDetElemId,fClusterNumber));
758 //_____________________________________________________________________________
760 AliMUONClusterFinderMLEM::ComputeCoefficients(AliMUONCluster& cluster,
764 /// Compute coefficients needed for MLEM algorithm
766 Int_t npadTot = cluster.Multiplicity();
767 Int_t nPix = fPixArray->GetLast()+1;
769 //memset(probi,0,nPix*sizeof(Double_t));
770 for (Int_t j = 0; j < npadTot*nPix; ++j) coef[j] = 0.;
771 for (Int_t j = 0; j < nPix; ++j) probi[j] = 0.;
773 Int_t mult = cluster.Multiplicity();
774 for ( Int_t j = 0; j < mult; ++j )
776 AliMUONPad* pad = cluster.Pad(j);
779 for ( Int_t ipix = 0; ipix < nPix; ++ipix )
781 Int_t indx1 = indx + ipix;
782 //if (pad->Status() < 0)
783 if (pad->Status() != fgkZero)
788 AliMUONPad* pixPtr = Pixel(ipix);
789 // coef is the charge (given by Mathieson integral) on pad, assuming
790 // the Mathieson is center at pixel.
791 coef[indx1] = fSplitter->ChargeIntegration(pixPtr->Coord(0), pixPtr->Coord(1), *pad);
792 // AliDebug(2,Form("pad=(%d,%d,%e,%e,%e,%e) pix=(%e,%e,%e,%e) coef %e",
793 // pad->Ix(),pad->Iy(),
794 // pad->X(),pad->Y(),
795 // pad->DX(),pad->DY(),
796 // pixPtr->Coord(0),pixPtr->Coord(1),
797 // pixPtr->Size(0),pixPtr->Size(1),
800 probi[ipix] += coef[indx1];
805 //_____________________________________________________________________________
806 Bool_t AliMUONClusterFinderMLEM::MainLoop(AliMUONCluster& cluster, Int_t iSimple)
808 /// Repeat MLEM algorithm until pixel size becomes sufficiently small
810 // AliCodeTimerAuto("")
812 Int_t nPix = fPixArray->GetLast()+1;
814 AliDebug(2,Form("nPix=%d iSimple=%d, precluster=",nPix,iSimple));
815 //StdoutToAliDebug(2,cluster.Print("full"););
819 AliDebug(1,"No pixels, why am I here then ?");
822 AddVirtualPad(cluster); // add virtual pads if necessary
824 Int_t npadTot = cluster.Multiplicity();
826 for (Int_t i = 0; i < npadTot; ++i)
828 //if (cluster.Pad(i)->Status() == 0) ++npadOK;
829 if (cluster.Pad(i)->Status() == fgkZero) ++npadOK;
833 Double_t* probi(0x0);
834 Int_t lc(0); // loop counter
836 //Plot("mlem.start");
837 AliMUONPad* pixPtr = Pixel(0);
838 Double_t xylim[4] = {pixPtr->X(), -pixPtr->X(), pixPtr->Y(), -pixPtr->Y()};
845 AliDebug(2,Form("lc %d nPix %d(%d) npadTot %d npadOK %d",lc,nPix,fPixArray->GetLast()+1,npadTot,npadOK));
846 AliDebug(2,Form("EVT%d PixArray=",fEventNumber));
847 //StdoutToAliDebug(2,fPixArray->Print("","full"));
849 coef = new Double_t [npadTot*nPix];
850 probi = new Double_t [nPix];
852 // Calculate coefficients and pixel visibilities
853 ComputeCoefficients(cluster,coef,probi);
855 for (Int_t ipix = 0; ipix < nPix; ++ipix)
857 if (probi[ipix] < 0.01)
859 AliMUONPad* pixel = Pixel(ipix);
860 AliDebug(2,Form("Setting the following pixel to invisible as its probi<0.01:"));
861 //StdoutToAliDebug(2,cout << Form(" -- ipix %3d --- "); pixel->Print(););
862 pixel->SetCharge(0); // "invisible" pixel
867 Mlem(cluster,coef, probi, 15);
869 // Find histogram limits for the 1'st pass only - for others computed below
871 for ( Int_t ipix = 1; ipix < nPix; ++ipix )
873 pixPtr = Pixel(ipix);
874 for ( Int_t i = 0; i < 2; ++i )
877 if (pixPtr->Coord(i) < xylim[indx]) xylim[indx] = pixPtr->Coord(i);
878 else if (-pixPtr->Coord(i) < xylim[indx+1]) xylim[indx+1] = -pixPtr->Coord(i);
881 } else pixPtr = Pixel(0);
883 for (Int_t i = 0; i < 4; i++)
885 xylim[i] -= pixPtr->Size(i/2);
888 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
889 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
891 //StdoutToAliDebug(2,cout << "pixel used for nx,ny computation : "; pixPtr->Print(););
892 AliDebug(2,Form("lc %d pixPtr size = %e,%e nx,ny=%d,%d xylim=%e,%e,%e,%e",
893 lc,pixPtr->Size(0),pixPtr->Size(1),nx,ny,
894 xylim[0],-xylim[1],xylim[2],-xylim[3]
897 fHistMlem = new TH2D("mlem","mlem",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
899 for (Int_t ipix = 0; ipix < nPix; ++ipix)
901 AliMUONPad* pixPtr2 = Pixel(ipix);
902 fHistMlem->Fill(pixPtr2->Coord(0),pixPtr2->Coord(1),pixPtr2->Charge());
905 // Check if the total charge of pixels is too low
907 for ( Int_t i = 0; i < nPix; ++i)
909 qTot += Pixel(i)->Charge();
912 if ( qTot < 1.e-4 || ( npadOK < 3 && qTot < 7 ) )
914 AliDebug(1,Form("Deleting the above cluster (charge %e too low, npadOK=%d)",qTot,npadOK));
918 for ( Int_t i = 0; i < npadTot; ++i)
920 AliMUONPad* pad = cluster.Pad(i);
921 //if ( pad->Status() == 0) pad->SetStatus(-1);
922 if ( pad->Status() == fgkZero) pad->SetStatus(fgkOver);
929 // Simple cluster - skip further passes thru EM-procedure
937 // Calculate position of the center-of-gravity around the maximum pixel
941 if (TMath::Min(pixPtr->Size(0),pixPtr->Size(1)) < 0.07 &&
942 pixPtr->Size(0) > pixPtr->Size(1)) break;
944 // Sort pixels according to the charge
945 MaskPeaks(1); // mask local maxima
947 MaskPeaks(0); // unmask local maxima
948 Double_t pixMin = 0.01*Pixel(0)->Charge();
949 pixMin = TMath::Min(pixMin,50.);
951 // Decrease pixel size and shift pixels to make them centered at
953 Int_t indx = (pixPtr->Size(0)>pixPtr->Size(1)) ? 0 : 1;
956 Double_t shift[2] = { 0.0, 0.0 };
957 for (Int_t i = 0; i < 4; ++i) xylim[i] = 999;
960 for (Int_t ipix = 0; ipix < nPix1; ++ipix)
962 AliMUONPad* pixPtr2 = Pixel(ipix);
963 if ( nPix >= npadOK // too many pixels already
965 ((pixPtr2->Charge() < pixMin) && (pixPtr2->Status() != fgkMustKeep)) // too low charge
971 for (Int_t i = 0; i < 2; ++i)
975 pixPtr2->SetCharge(10);
976 pixPtr2->SetSize(indx, pixPtr2->Size(indx)/2);
977 width = -pixPtr2->Size(indx);
978 pixPtr2->Shift(indx, width);
979 // Shift pixel position
983 for (Int_t j = 0; j < 2; ++j)
985 shift[j] = pixPtr2->Coord(j) - xyCOG[j];
986 shift[j] -= ((Int_t)(shift[j]/pixPtr2->Size(j)/2))*pixPtr2->Size(j)*2;
989 pixPtr2->Shift(0, -shift[0]);
990 pixPtr2->Shift(1, -shift[1]);
993 else if (nPix < npadOK)
995 pixPtr2 = new AliMUONPad(*pixPtr2);
996 pixPtr2->Shift(indx, -2*width);
997 pixPtr2->SetStatus(fgkZero);
998 fPixArray->Add(pixPtr2);
1001 else continue; // skip adjustment of histo limits
1002 for (Int_t j = 0; j < 4; ++j)
1004 xylim[j] = TMath::Min (xylim[j], (j%2 ? -1 : 1)*pixPtr2->Coord(j/2));
1006 } // for (Int_t i=0; i<2;
1007 } // for (Int_t ipix=0;
1009 fPixArray->Compress();
1011 AliDebug(2,Form("After shift:"));
1012 //StdoutToAliDebug(2,fPixArray->Print("","full"););
1013 //Plot(Form("mlem.lc%d",lc+1));
1015 AliDebug(2,Form(" xyCOG=%9.6f %9.6f xylim=%9.6f,%9.6f,%9.6f,%9.6f",
1018 xylim[2],xylim[3]));
1022 AliMUONPad* pixPtr2 = Pixel(0);
1023 // add pixels if the maximum is at the limit of pixel area:
1024 // start from Y-direction
1026 for (Int_t i = 3; i > -1; --i)
1028 if (nPix < npadOK &&
1029 TMath::Abs((i%2 ? -1 : 1)*xylim[i]-xyCOG[i/2]) < pixPtr2->Size(i/2))
1031 //AliMUONPad* p = static_cast<AliMUONPad*>(pixPtr->Clone());
1032 AliMUONPad* p = new AliMUONPad(*pixPtr2);
1033 p->SetCoord(i/2, xyCOG[i/2]+(i%2 ? 2:-2)*pixPtr2->Size(i/2));
1034 xylim[i] = p->Coord(i/2) * (i%2 ? -1 : 1); // update histo limits
1035 j = TMath::Even (i/2);
1036 p->SetCoord(j, xyCOG[j]);
1037 AliDebug(2,Form("Adding pixel on the edge (i=%d) ",i));
1038 //StdoutToAliDebug(2,cout << " ---- ";
1039 // p->Print("corners"););
1045 nPix = fPixArray->GetEntriesFast();
1050 AliDebug(2,Form("At the end of while loop nPix=%d : ",fPixArray->GetLast()+1));
1051 //StdoutToAliDebug(2,fPixArray->Print("","full"););
1053 // remove pixels with low signal or low visibility
1054 // Cuts are empirical !!!
1055 Double_t thresh = TMath::Max (fHistMlem->GetMaximum()/100.,1.);
1056 thresh = TMath::Min (thresh,50.);
1057 Double_t charge = 0;
1059 // Mark pixels which should be removed
1060 for (Int_t i = 0; i < nPix; ++i)
1062 AliMUONPad* pixPtr2 = Pixel(i);
1063 charge = pixPtr2->Charge();
1064 if (charge < thresh)
1066 pixPtr2->SetCharge(-charge);
1070 // Move charge of removed pixels to their nearest neighbour (to keep total charge the same)
1072 for (Int_t i = 0; i < nPix; ++i)
1074 AliMUONPad* pixPtr2 = Pixel(i);
1075 charge = pixPtr2->Charge();
1076 if (charge > 0) continue;
1077 near = FindNearest(pixPtr2);
1078 pixPtr2->SetCharge(0);
1079 probi[i] = 0; // make it "invisible"
1080 AliMUONPad* pnear = Pixel(near);
1081 pnear->SetCharge(pnear->Charge() + (-charge));
1083 Mlem(cluster,coef,probi,2);
1085 AliDebug(2,Form("Before splitting nPix=%d EVT %d DE %d",fPixArray->GetLast()+1,fEventNumber,fDetElemId));
1086 //StdoutToAliDebug(2,fPixArray->Print("","full"););
1087 //Plot("mlem.beforesplit");
1090 for (Int_t i = 0; i < nPix; ++i)
1092 AliMUONPad* pixPtr2 = Pixel(i);
1093 Int_t ix = fHistMlem->GetXaxis()->FindBin(pixPtr2->Coord(0));
1094 Int_t iy = fHistMlem->GetYaxis()->FindBin(pixPtr2->Coord(1));
1095 fHistMlem->SetBinContent(ix, iy, pixPtr2->Charge());
1098 // Try to split into clusters
1100 if (fHistMlem->GetSum() < 1)
1106 fSplitter->Split(cluster,fHistMlem,coef,fClusterList);
1111 fPixArray->Delete();
1116 //_____________________________________________________________________________
1117 void AliMUONClusterFinderMLEM::MaskPeaks(Int_t mask)
1119 /// Mask/unmask pixels corresponding to local maxima (add/subtract 10000 to their charge
1120 /// - to avoid loosing low charge pixels after sorting)
1122 for (Int_t i = 0; i < fPixArray->GetEntriesFast(); ++i) {
1123 AliMUONPad* pix = Pixel(i);
1124 if (pix->Status() == fgkMustKeep) {
1125 if (mask == 1) pix->SetCharge(pix->Charge()+10000.);
1126 else pix->SetCharge(pix->Charge()-10000.);
1131 //_____________________________________________________________________________
1132 void AliMUONClusterFinderMLEM::Mlem(AliMUONCluster& cluster,
1133 const Double_t* coef, Double_t* probi,
1136 /// Use MLEM to find pixel charges
1138 Int_t nPix = fPixArray->GetEntriesFast();
1140 Int_t npad = cluster.Multiplicity();
1142 Double_t* probi1 = new Double_t[nPix];
1143 Double_t probMax = TMath::MaxElement(nPix,probi);
1145 for (Int_t iter = 0; iter < nIter; ++iter)
1148 for (Int_t ipix = 0; ipix < nPix; ++ipix)
1150 Pixel(ipix)->SetChargeBackup(0);
1151 // Correct each pixel
1153 if (probi[ipix] < 0.01) continue; // skip "invisible" pixel
1155 probi1[ipix] = probMax;
1156 for (Int_t j = 0; j < npad; ++j)
1158 AliMUONPad* pad = cluster.Pad(j);
1159 //if (pad->Status() < 0) continue;
1160 if (pad->Status() != fgkZero) continue;
1162 Int_t indx1 = j*nPix;
1163 Int_t indx = indx1 + ipix;
1164 // Calculate expectation
1165 for (Int_t i = 0; i < nPix; ++i)
1167 sum1 += Pixel(i)->Charge()*coef[indx1+i];
1168 //cout << i << " " << Pixel(i)->Charge() << " " << coef[indx1+i] << endl;
1170 if ( pad->IsSaturated() && sum1 > pad->Charge() )
1172 // correct for pad charge overflows
1173 probi1[ipix] -= coef[indx];
1175 //sum1 = pad->Charge();
1178 if (sum1 > 1.e-6) sum += pad->Charge()*coef[indx]/sum1;
1179 } // for (Int_t j=0;
1180 AliMUONPad* pixPtr = Pixel(ipix);
1181 if (probi1[ipix] > 1.e-6)
1183 //AZ pixPtr->SetCharge(pixPtr->Charge()*sum/probi1[ipix]);
1184 pixPtr->SetChargeBackup(pixPtr->Charge()*sum/probi1[ipix]);
1186 //cout << " xxx " << ipix << " " << pixPtr->Charge() << " " << pixPtr->ChargeBackup() << " " << sum << " " << probi1[ipix] << endl;
1187 } // for (Int_t ipix=0;
1189 for (Int_t i = 0; i < nPix; ++i) {
1190 AliMUONPad* pixPtr = Pixel(i);
1191 pixPtr->RevertCharge();
1192 qTot += pixPtr->Charge();
1195 // Can happen in clusters with large number of overflows - speeding up
1199 } // for (Int_t iter=0;
1203 //_____________________________________________________________________________
1204 void AliMUONClusterFinderMLEM::FindCOG(Double_t *xyc)
1206 /// Calculate position of the center-of-gravity around the maximum pixel
1208 Int_t ixmax, iymax, ix, nsumx=0, nsumy=0, nsum=0;
1209 Int_t i1 = -9, j1 = -9;
1210 fHistMlem->GetMaximumBin(ixmax,iymax,ix);
1211 Int_t nx = fHistMlem->GetNbinsX();
1212 Int_t ny = fHistMlem->GetNbinsY();
1213 Double_t thresh = fHistMlem->GetMaximum()/10;
1214 Double_t x, y, cont, xq=0, yq=0, qq=0;
1216 Int_t ie = TMath::Min(ny,iymax+1), je = TMath::Min(nx,ixmax+1);
1217 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
1218 y = fHistMlem->GetYaxis()->GetBinCenter(i);
1219 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
1220 cont = fHistMlem->GetCellContent(j,i);
1221 if (cont < thresh) continue;
1222 if (i != i1) {i1 = i; nsumy++;}
1223 if (j != j1) {j1 = j; nsumx++;}
1224 x = fHistMlem->GetXaxis()->GetBinCenter(j);
1233 Int_t i2 = 0, j2 = 0;
1236 // one bin in Y - add one more (with the largest signal)
1237 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
1238 if (i == iymax) continue;
1239 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
1240 cont = fHistMlem->GetCellContent(j,i);
1243 x = fHistMlem->GetXaxis()->GetBinCenter(j);
1244 y = fHistMlem->GetYaxis()->GetBinCenter(i);
1253 if (i2 != i1) nsumy++;
1254 if (j2 != j1) nsumx++;
1256 } // if (nsumy == 1)
1259 // one bin in X - add one more (with the largest signal)
1261 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
1262 if (j == ixmax) continue;
1263 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
1264 cont = fHistMlem->GetCellContent(j,i);
1267 x = fHistMlem->GetXaxis()->GetBinCenter(j);
1268 y = fHistMlem->GetYaxis()->GetBinCenter(i);
1277 if (i2 != i1) nsumy++;
1278 if (j2 != j1) nsumx++;
1280 } // if (nsumx == 1)
1282 xyc[0] = xq/qq; xyc[1] = yq/qq;
1283 AliDebug(2,Form("x,y COG = %e,%e",xyc[0],xyc[1]));
1286 //_____________________________________________________________________________
1287 Int_t AliMUONClusterFinderMLEM::FindNearest(const AliMUONPad *pixPtr0)
1289 /// Find the pixel nearest to the given one
1290 /// (algorithm may be not very efficient)
1292 Int_t nPix = fPixArray->GetEntriesFast(), imin = 0;
1293 Double_t rmin = 99999, dx = 0, dy = 0, r = 0;
1294 Double_t xc = pixPtr0->Coord(0), yc = pixPtr0->Coord(1);
1297 for (Int_t i = 0; i < nPix; ++i) {
1298 pixPtr = (AliMUONPad*) fPixArray->UncheckedAt(i);
1299 if (pixPtr == pixPtr0 || pixPtr->Charge() < 0.5) continue;
1300 dx = (xc - pixPtr->Coord(0)) / pixPtr->Size(0);
1301 dy = (yc - pixPtr->Coord(1)) / pixPtr->Size(1);
1302 r = dx *dx + dy * dy;
1303 if (r < rmin) { rmin = r; imin = i; }
1308 //_____________________________________________________________________________
1310 AliMUONClusterFinderMLEM::Paint(Option_t*)
1312 /// Paint cluster and pixels
1314 AliMpArea area(fPreCluster->Area());
1316 gPad->Range(area.LeftBorder(),area.DownBorder(),area.RightBorder(),area.UpBorder());
1318 gVirtualX->SetFillStyle(1001);
1319 gVirtualX->SetFillColor(3);
1320 gVirtualX->SetLineColor(3);
1324 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1326 AliMUONPad* pixel = Pixel(i);
1328 gPad->PaintBox(pixel->Coord(0)-pixel->Size(0)*s,
1329 pixel->Coord(1)-pixel->Size(1)*s,
1330 pixel->Coord(0)+pixel->Size(0)*s,
1331 pixel->Coord(1)+pixel->Size(1)*s);
1333 // for ( Int_t sign = -1; sign < 2; sign +=2 )
1335 // gPad->PaintLine(pixel->Coord(0) - pixel->Size(0),
1336 // pixel->Coord(1) + sign*pixel->Size(1),
1337 // pixel->Coord(0) + pixel->Size(0),
1338 // pixel->Coord(1) - sign*pixel->Size(1)
1344 gVirtualX->SetFillStyle(0);
1346 fPreCluster->Paint();
1348 gVirtualX->SetLineColor(1);
1349 gVirtualX->SetLineWidth(2);
1350 gVirtualX->SetFillStyle(0);
1351 gVirtualX->SetTextColor(1);
1352 gVirtualX->SetTextAlign(22);
1354 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1356 AliMUONPad* pixel = Pixel(i);
1357 gPad->PaintBox(pixel->Coord(0)-pixel->Size(0),
1358 pixel->Coord(1)-pixel->Size(1),
1359 pixel->Coord(0)+pixel->Size(0),
1360 pixel->Coord(1)+pixel->Size(1));
1361 gVirtualX->SetTextSize(pixel->Size(0)*60);
1363 gPad->PaintText(pixel->Coord(0),pixel->Coord(1),Form("%d",(Int_t)(pixel->Charge())));
1367 //_____________________________________________________________________________
1368 Int_t AliMUONClusterFinderMLEM::FindLocalMaxima(TObjArray *pixArray, Int_t *localMax, Double_t *maxVal)
1370 /// Find local maxima in pixel space for large preclusters in order to
1371 /// try to split them into smaller pieces (to speed up the MLEM procedure)
1372 /// or to find additional fitting seeds if clusters were not completely resolved
1374 AliDebug(1,Form("nPix=%d",pixArray->GetLast()+1));
1376 Double_t xylim[4] = {999, 999, 999, 999};
1378 Int_t nPix = pixArray->GetEntriesFast();
1379 AliMUONPad *pixPtr = 0;
1380 for (Int_t ipix = 0; ipix < nPix; ++ipix) {
1381 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
1382 for (Int_t i = 0; i < 4; ++i)
1383 xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
1385 for (Int_t i = 0; i < 4; ++i) xylim[i] -= pixPtr->Size(i/2);
1387 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
1388 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
1389 if (pixArray == fPixArray) fHistAnode = new TH2D("anode","anode",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1390 else fHistAnode = new TH2D("anode1","anode1",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1391 for (Int_t ipix = 0; ipix < nPix; ++ipix) {
1392 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
1393 fHistAnode->Fill(pixPtr->Coord(0), pixPtr->Coord(1), pixPtr->Charge());
1395 // if (fDraw && pixArray == fPixArray) fDraw->DrawHist("c2", hist);
1397 Int_t nMax = 0, indx, nxy = ny * nx;
1398 Int_t *isLocalMax = new Int_t[nxy];
1399 for (Int_t i = 0; i < nxy; ++i) isLocalMax[i] = 0;
1401 for (Int_t i = 1; i <= ny; ++i) {
1403 for (Int_t j = 1; j <= nx; ++j) {
1404 if (fHistAnode->GetCellContent(j,i) < 0.5) continue;
1405 //if (isLocalMax[indx+j-1] < 0) continue;
1406 if (isLocalMax[indx+j-1] != 0) continue;
1407 FlagLocalMax(fHistAnode, i, j, isLocalMax);
1411 for (Int_t i = 1; i <= ny; ++i) {
1413 for (Int_t j = 1; j <= nx; ++j) {
1414 if (isLocalMax[indx+j-1] > 0) {
1415 localMax[nMax] = indx + j - 1;
1416 maxVal[nMax++] = fHistAnode->GetCellContent(j,i);
1417 ((AliMUONPad*)fSplitter->BinToPix(fHistAnode, j, i))->SetStatus(fgkMustKeep);
1418 if (nMax > 99) break;
1422 AliError(" Too many local maxima !!!");
1426 if (fDebug) cout << " Local max: " << nMax << endl;
1427 delete [] isLocalMax;
1431 //_____________________________________________________________________________
1432 void AliMUONClusterFinderMLEM::FlagLocalMax(TH2D *hist, Int_t i, Int_t j, Int_t *isLocalMax)
1434 /// Flag pixels (whether or not local maxima)
1436 Int_t nx = hist->GetNbinsX();
1437 Int_t ny = hist->GetNbinsY();
1438 Int_t cont = TMath::Nint (hist->GetCellContent(j,i));
1439 Int_t cont1 = 0, indx = (i-1)*nx+j-1, indx1 = 0, indx2 = 0;
1441 Int_t ie = i + 2, je = j + 2;
1442 for (Int_t i1 = i-1; i1 < ie; ++i1) {
1443 if (i1 < 1 || i1 > ny) continue;
1444 indx1 = (i1 - 1) * nx;
1445 for (Int_t j1 = j-1; j1 < je; ++j1) {
1446 if (j1 < 1 || j1 > nx) continue;
1447 if (i == i1 && j == j1) continue;
1448 indx2 = indx1 + j1 - 1;
1449 cont1 = TMath::Nint (hist->GetCellContent(j1,i1));
1450 if (cont < cont1) { isLocalMax[indx] = -1; return; }
1451 else if (cont > cont1) isLocalMax[indx2] = -1;
1452 else { // the same charge
1453 isLocalMax[indx] = 1;
1454 if (isLocalMax[indx2] == 0) {
1455 FlagLocalMax(hist, i1, j1, isLocalMax);
1456 if (isLocalMax[indx2] < 0) { isLocalMax[indx] = -1; return; }
1457 else isLocalMax[indx2] = -1;
1462 isLocalMax[indx] = 1; // local maximum
1465 //_____________________________________________________________________________
1466 void AliMUONClusterFinderMLEM::FindCluster(AliMUONCluster& cluster,
1467 const Int_t *localMax, Int_t iMax)
1469 /// Find pixel cluster around local maximum \a iMax and pick up pads
1470 /// overlapping with it
1473 TCanvas* c = new TCanvas("Anode","Anode",800,600);
1475 hist->Draw("lego1Fb"); // debug
1480 Int_t nx = fHistAnode->GetNbinsX();
1481 Int_t ny = fHistAnode->GetNbinsY();
1482 Int_t ic = localMax[iMax] / nx + 1;
1483 Int_t jc = localMax[iMax] % nx + 1;
1484 Int_t nxy = ny * nx;
1485 Bool_t *used = new Bool_t[nxy];
1486 for (Int_t i = 0; i < nxy; ++i) used[i] = kFALSE;
1488 // Drop all pixels from the array - pick up only the ones from the cluster
1489 fPixArray->Delete();
1491 Double_t wx = fHistAnode->GetXaxis()->GetBinWidth(1)/2;
1492 Double_t wy = fHistAnode->GetYaxis()->GetBinWidth(1)/2;
1493 Double_t yc = fHistAnode->GetYaxis()->GetBinCenter(ic);
1494 Double_t xc = fHistAnode->GetXaxis()->GetBinCenter(jc);
1495 Double_t cont = fHistAnode->GetCellContent(jc,ic);
1496 fPixArray->Add(new AliMUONPad (xc, yc, wx, wy, cont));
1497 used[(ic-1)*nx+jc-1] = kTRUE;
1498 AddBinSimple(fHistAnode, ic, jc);
1499 //fSplitter->AddBin(hist, ic, jc, 1, used, (TObjArray*)0); // recursive call
1501 Int_t nPix = fPixArray->GetEntriesFast();
1502 Int_t npad = cluster.Multiplicity();
1504 for (Int_t i = 0; i < nPix; ++i)
1506 AliMUONPad* pixPtr = Pixel(i);
1507 pixPtr->SetSize(0,wx);
1508 pixPtr->SetSize(1,wy);
1511 // Pick up pads which overlap with found pixels
1512 for (Int_t i = 0; i < npad; ++i)
1514 //cluster.Pad(i)->SetStatus(-1);
1515 cluster.Pad(i)->SetStatus(fgkOver); // just the dirty trick
1518 for (Int_t i = 0; i < nPix; ++i)
1520 AliMUONPad* pixPtr = Pixel(i);
1521 for (Int_t j = 0; j < npad; ++j)
1523 AliMUONPad* pad = cluster.Pad(j);
1524 //if (pad->Status() == 0) continue;
1525 if (pad->Status() == fgkZero) continue;
1526 if ( Overlap(*pad,*pixPtr) )
1528 //pad->SetStatus(0);
1529 pad->SetStatus(fgkZero);
1530 if (fDebug) { cout << j << " "; pad->Print("full"); }
1538 //_____________________________________________________________________________
1540 AliMUONClusterFinderMLEM::AddBinSimple(TH2D *hist, Int_t ic, Int_t jc)
1542 /// Add adjacent bins (+-1 in X and Y) to the cluster
1544 Int_t nx = hist->GetNbinsX();
1545 Int_t ny = hist->GetNbinsY();
1546 Double_t cont1, cont = hist->GetCellContent(jc,ic);
1547 AliMUONPad *pixPtr = 0;
1549 Int_t ie = TMath::Min(ic+1,ny), je = TMath::Min(jc+1,nx);
1550 for (Int_t i = TMath::Max(ic-1,1); i <= ie; ++i) {
1551 for (Int_t j = TMath::Max(jc-1,1); j <= je; ++j) {
1552 cont1 = hist->GetCellContent(j,i);
1553 if (cont1 > cont) continue;
1554 if (cont1 < 0.5) continue;
1555 pixPtr = new AliMUONPad (hist->GetXaxis()->GetBinCenter(j),
1556 hist->GetYaxis()->GetBinCenter(i), 0, 0, cont1);
1557 fPixArray->Add(pixPtr);
1562 //_____________________________________________________________________________
1563 AliMUONClusterFinderMLEM&
1564 AliMUONClusterFinderMLEM::operator=(const AliMUONClusterFinderMLEM& rhs)
1566 /// Protected assignement operator
1568 if (this == &rhs) return *this;
1570 AliFatal("Not implemented.");
1575 //_____________________________________________________________________________
1576 void AliMUONClusterFinderMLEM::AddVirtualPad(AliMUONCluster& cluster)
1578 /// Add virtual pad (with small charge) to improve fit for clusters
1579 /// with number of pads == 2 per direction
1581 // Find out non-bending and bending planes
1582 Int_t nonb[2] = {1, 0}; // non-bending and bending cathodes
1584 TVector2 dim0 = cluster.MinPadDimensions(0, 0, kTRUE);
1585 TVector2 dim1 = cluster.MinPadDimensions(1, 0, kTRUE);
1586 if (dim0.X() < dim1.X() - fgkDistancePrecision) {
1591 Bool_t same = kFALSE;
1592 if (TMath::Abs(dim0.Y()-dim1.Y()) < fgkDistancePrecision) same = kTRUE; // the same pad size on both planes
1595 Bool_t check[2] = {kFALSE, kFALSE};
1597 nxy[0] = nxy[1] = 0;
1598 for (Int_t inb = 0; inb < 2; ++inb) {
1599 cn = cluster.NofPads(nonb[inb], 0, kTRUE);
1600 if (inb == 0 && AliMp::PairFirst(cn) == 2) check[inb] = kTRUE; // check non-bending plane
1601 else if (inb == 1 && AliMp::PairSecond(cn) == 2) check[inb] = kTRUE; // check bending plane
1603 nxy[0] = TMath::Max (nxy[0], AliMp::PairFirst(cn));
1604 nxy[1] = TMath::Max (nxy[1], AliMp::PairSecond(cn));
1605 if (inb == 0 && nxy[0] < 2) nonb[inb] = !nonb[inb];
1606 else if (inb == 1 && AliMp::PairSecond(cn) < 2) nonb[inb] = !nonb[inb];
1610 if (nxy[0] > 2) check[0] = kFALSE;
1611 if (nxy[1] > 2) check[1] = kFALSE;
1613 if (!check[0] && !check[1]) return;
1615 for (Int_t inb = 0; inb < 2; ++inb) {
1616 if (!check[inb]) continue;
1618 // Find pads with maximum and next to maximum charges
1619 Int_t maxPads[2] = {-1, -1};
1620 Double_t amax[2] = {0};
1621 Int_t mult = cluster.Multiplicity();
1622 for (Int_t j = 0; j < mult; ++j) {
1623 AliMUONPad *pad = cluster.Pad(j);
1624 if (pad->Cathode() != nonb[inb]) continue;
1625 for (Int_t j2 = 0; j2 < 2; ++j2) {
1626 if (pad->Charge() > amax[j2]) {
1627 if (j2 == 0) { amax[1] = amax[0]; maxPads[1] = maxPads[0]; }
1628 amax[j2] = pad->Charge();
1635 // Find min and max dimensions of the cluster
1636 Double_t limits[2] = {9999, -9999};
1637 for (Int_t j = 0; j < mult; ++j) {
1638 AliMUONPad *pad = cluster.Pad(j);
1639 if (pad->Cathode() != nonb[inb]) continue;
1640 if (pad->Coord(inb) < limits[0]) limits[0] = pad->Coord(inb);
1641 if (pad->Coord(inb) > limits[1]) limits[1] = pad->Coord(inb);
1644 // Loop over max and next to max pads
1645 Bool_t add = kFALSE;
1647 for (Int_t j = 0; j < 2; ++j) {
1649 if (maxPads[j] < 0) continue;
1651 if (amax[1] / amax[0] < 0.5) break;
1653 // Check if pad at the cluster limit
1654 AliMUONPad *pad = cluster.Pad(maxPads[j]);
1656 if (TMath::Abs(pad->Coord(inb)-limits[0]) < fgkDistancePrecision) idir = -1;
1657 else if (TMath::Abs(pad->Coord(inb)-limits[1]) < fgkDistancePrecision) idir = 1;
1659 //cout << " *** Pad not at the cluster limit: " << j << endl;
1662 if (j == 1 && idir == idirAdd) break; // this pad has been already added
1664 // Add pad (if it exists)
1666 if (inb == 0) pos.Set (pad->X() + idir * (pad->DX()+fgkDistancePrecision), pad->Y());
1667 else pos.Set (pad->X(), pad->Y() + idir * (pad->DY()+fgkDistancePrecision));
1668 //AliMpPad mppad = fkSegmentation[nonb[inb]]->PadByPosition(pos,kTRUE);
1669 AliMpPad mppad = fkSegmentation[nonb[inb]]->PadByPosition(pos.X(), pos.Y(),kFALSE);
1670 if (!mppad.IsValid()) continue; // non-existing pad
1671 AliMUONPad muonPad(fDetElemId, nonb[inb], mppad.GetIx(), mppad.GetIy(),
1672 mppad.GetPositionX(), mppad.GetPositionY(),
1673 mppad.GetDimensionX(), mppad.GetDimensionY(), 0);
1674 if (inb == 0) muonPad.SetCharge(TMath::Min (amax[j]/100, 5.));
1675 //else muonPad.SetCharge(TMath::Min (amax[j]/15, fgkZeroSuppression));
1676 else muonPad.SetCharge(TMath::Min (amax[j]/15, 6.));
1677 if (muonPad.Charge() < 1.) muonPad.SetCharge(1.);
1678 muonPad.SetReal(kFALSE);
1679 if (fDebug) printf(" ***** Add virtual pad in %d direction ***** %f %f %f %3d %3d %f %f \n",
1680 inb, muonPad.Charge(), muonPad.X(), muonPad.Y(), muonPad.Ix(),
1681 muonPad.Iy(), muonPad.DX(), muonPad.DY());
1682 cluster.AddPad(muonPad); // add pad to the cluster
1689 //_____________________________________________________________________________
1690 void AliMUONClusterFinderMLEM::PadsInXandY(AliMUONCluster& cluster,
1691 Int_t &nInX, Int_t &nInY) const
1693 /// Find number of pads in X and Y-directions (excluding virtual ones and
1696 //Int_t statusToTest = 1;
1697 Int_t statusToTest = fgkUseForFit;
1699 //if ( nInX < 0 ) statusToTest = 0;
1700 if ( nInX < 0 ) statusToTest = fgkZero;
1702 Bool_t mustMatch(kTRUE);
1704 Long_t cn = cluster.NofPads(statusToTest,mustMatch);
1706 nInX = AliMp::PairFirst(cn);
1707 nInY = AliMp::PairSecond(cn);
1710 //_____________________________________________________________________________
1711 void AliMUONClusterFinderMLEM::RemovePixel(Int_t i)
1713 /// Remove pixel at index i
1714 AliMUONPad* pixPtr = Pixel(i);
1715 fPixArray->RemoveAt(i);
1719 //_____________________________________________________________________________
1720 void AliMUONClusterFinderMLEM::Simple(AliMUONCluster& cluster)
1722 /// Process simple cluster (small number of pads) without EM-procedure
1724 Int_t nForFit = 1, clustFit[1] = {0}, nfit;
1725 Double_t parOk[3] = {0.};
1726 TObjArray *clusters[1];
1727 clusters[0] = fPixArray;
1729 AliDebug(1,Form("nPix=%d",fPixArray->GetLast()+1));
1731 Int_t mult = cluster.Multiplicity();
1732 for (Int_t i = 0; i < mult; ++i)
1734 AliMUONPad* pad = cluster.Pad(i);
1736 if ( pad->IsSaturated())
1745 if (!pad->IsSaturated()) pad->SetStatus(fgkUseForFit);
1747 nfit = fSplitter->Fit(cluster,1, nForFit, clustFit, clusters, parOk, fClusterList, fHistMlem);
1750 //_____________________________________________________________________________
1752 AliMUONClusterFinderMLEM::Pixel(Int_t i) const
1754 /// Returns pixel at index i
1755 return static_cast<AliMUONPad*>(fPixArray->UncheckedAt(i));
1758 //_____________________________________________________________________________
1760 AliMUONClusterFinderMLEM::Print(Option_t* what) const
1763 TString swhat(what);
1765 if ( swhat.Contains("precluster") )
1767 if ( fPreCluster) fPreCluster->Print();