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 /// \class AliMUONClusterFinderMLEM
20 /// Clusterizer class based on the Expectation-Maximization algorithm
22 /// Pre-clustering is handled by AliMUONPreClusterFinder
23 /// From a precluster a pixel array is built, and from this pixel array
24 /// a list of clusters is output (using AliMUONClusterSplitterMLEM).
26 /// \author Laurent Aphecetche (for the "new" C++ structure) and
27 /// Alexander Zinchenko, JINR Dubna, for the hardcore of it ;-)
29 #include "AliMUONClusterFinderMLEM.h"
31 #include "AliMUONCluster.h"
32 #include "AliMUONClusterSplitterMLEM.h"
33 #include "AliMUONDigit.h"
34 #include "AliMUONPad.h"
35 #include "AliMUONPreClusterFinder.h"
37 #include "AliMpVPadIterator.h"
38 #include "AliMpVSegmentation.h"
39 #include "AliRunLoader.h"
41 #include <Riostream.h>
45 #include <TStopwatch.h>
50 ClassImp(AliMUONClusterFinderMLEM)
53 const Double_t AliMUONClusterFinderMLEM::fgkZeroSuppression = 6; // average zero suppression value
54 const Double_t AliMUONClusterFinderMLEM::fgkSaturation = 3000; // average saturation level
55 const Double_t AliMUONClusterFinderMLEM::fgkDistancePrecision = 1e-6; // (cm) used to check overlaps and so on
56 const TVector2 AliMUONClusterFinderMLEM::fgkIncreaseSize(-AliMUONClusterFinderMLEM::fgkDistancePrecision,-AliMUONClusterFinderMLEM::fgkDistancePrecision);
57 const TVector2 AliMUONClusterFinderMLEM::fgkDecreaseSize(AliMUONClusterFinderMLEM::fgkDistancePrecision,AliMUONClusterFinderMLEM::fgkDistancePrecision);
59 TMinuit* AliMUONClusterFinderMLEM::fgMinuit = 0x0;
61 //_____________________________________________________________________________
62 AliMUONClusterFinderMLEM::AliMUONClusterFinderMLEM(Bool_t plot)
63 : AliMUONVClusterFinder(),
64 fPreClusterFinder(new AliMUONPreClusterFinder),
73 fPixArray(new TObjArray(20)),
76 fTimers(new TObjArray(kLast)),
83 fSegmentation[1] = fSegmentation[0] = 0x0;
85 fPadBeg[0] = fPadBeg[1] = fCathBeg;
87 if (!fgMinuit) fgMinuit = new TMinuit(8);
89 fTimers->SetOwner(kTRUE);
91 for ( Int_t i = 0; i < kLast; ++i )
93 TStopwatch* t = new TStopwatch;
94 fTimers->AddLast(new TStopwatch);
100 //_____________________________________________________________________________
101 AliMUONClusterFinderMLEM::~AliMUONClusterFinderMLEM()
104 delete fgMinuit; fgMinuit = 0; delete fPixArray; fPixArray = 0;
106 delete fPreClusterFinder;
107 for ( Int_t i = 0; i < kLast; ++i )
109 AliInfo(Form("Timer %d",i));
114 AliInfo(Form("Total clusters %d AddVirtualPad needed %d",
115 fNClusters,fNAddVirtualPads));
118 //_____________________________________________________________________________
120 AliMUONClusterFinderMLEM::Prepare(const AliMpVSegmentation* segmentations[2],
121 TClonesArray* digits[2])
123 /// Prepare for clustering
125 for ( Int_t i = 0; i < 2; ++i )
127 fSegmentation[i] = segmentations[i];
130 // Find out the DetElemId
133 for ( Int_t i = 0; i < 2; ++i )
135 AliMUONDigit* d = static_cast<AliMUONDigit*>(digits[i]->First());
138 fDetElemId = d->DetElemId();
143 if ( fDetElemId < 0 )
145 AliWarning("Could not find DE. Probably no digits at all ?");
150 fSplitter = new AliMUONClusterSplitterMLEM(fDetElemId,fPixArray);
152 // find out current event number, and reset the cluster number
153 fEventNumber = AliRunLoader::GetRunLoader()->GetEventNumber();
155 fClusterList.Delete();
157 // AliDebug(3,Form("EVT %d DE %d",fEventNumber,fDetElemId));
159 return fPreClusterFinder->Prepare(segmentations,digits);
162 //_____________________________________________________________________________
164 AliMUONClusterFinderMLEM::NextCluster()
166 /// Return next cluster
170 // if the list of clusters is not void, pick one from there
171 if ( fClusterList.GetLast() >= 0 )
173 TObject* o = fClusterList.At(0);
174 fClusterList.RemoveAt(0);
175 return static_cast<AliMUONCluster*>(o);
178 //FIXME : at this point, must check whether we've used all the digits
179 //from precluster : if not, let the preclustering know about those unused
180 //digits, so it can reuse them
182 // if the cluster list is exhausted, we need to go to the next
183 // pre-cluster and treat it
185 fPreCluster = fPreClusterFinder->NextCluster();
193 fClusterList.Delete(); // reset the list of clusters for this pre-cluster
197 // WorkOnPreCluster may have used only part of the pads, so we check that
198 // now, and let the unused pads be reused by the preclustering...
200 for ( Int_t i = 0; i < fPreCluster->Multiplicity(); ++i )
202 AliMUONPad* pad = fPreCluster->Pad(i);
203 if ( !pad->IsUsed() )
205 fPreClusterFinder->UsePad(*pad);
209 return NextCluster();
212 //_____________________________________________________________________________
214 AliMUONClusterFinderMLEM::WorkOnPreCluster()
216 /// Starting from a precluster, builds a pixel array, and then
217 /// extract clusters from this array
219 AliMUONCluster* cluster = CheckPrecluster(*fPreCluster);
221 if (!cluster) return kFALSE;
223 BuildPixArray(*cluster);
225 if ( fPixArray->GetLast() < 0 )
227 AliDebug(1,"No pixel for the above cluster");
232 // Use MLEM for cluster finder
233 Int_t nMax = 1, localMax[100], maxPos[100];
234 Double_t maxVal[100];
236 if (cluster->Multiplicity() > 50)
238 nMax = FindLocalMaxima(fPixArray, localMax, maxVal);
243 TMath::Sort(nMax, maxVal, maxPos, kTRUE); // in decreasing order
246 Int_t iSimple = 0, nInX = -1, nInY;
248 PadsInXandY(*cluster,nInX, nInY);
250 if (nMax == 1 && nInX < 4 && nInY < 4)
252 iSimple = 1; //1; // simple cluster
255 for (Int_t i=0; i<nMax; ++i)
259 FindCluster(*cluster,localMax, maxPos[i]);
261 Timer(kMainLoop)->Start(kFALSE);
262 MainLoop(*cluster,iSimple);
263 Timer(kMainLoop)->Stop();
266 for (Int_t j=0; j<cluster->Multiplicity(); ++j)
268 AliMUONPad* pad = cluster->Pad(j);
269 if ( pad->Status() == 0 ) continue; // pad charge was not modified
271 pad->RevertCharge(); // use backup charge value
274 } // for (Int_t i=0; i<nMax;
275 if (nMax > 1) ((TH2D*) gROOT->FindObject("anode"))->Delete();
276 TH2D *mlem = (TH2D*) gROOT->FindObject("mlem");
277 if (mlem) mlem->Delete();
282 //_____________________________________________________________________________
284 AliMUONClusterFinderMLEM::Overlap(const AliMUONPad& pad, const AliMUONPad& pixel)
286 /// Check if the pad and the pixel overlaps
288 // make a fake pad from the pixel
289 AliMUONPad tmp(pad.DetElemId(),pad.Cathode(),pad.Ix(),pad.Iy(),
290 pixel.Coord(0),pixel.Coord(1),
291 pixel.Size(0),pixel.Size(1),0);
293 return AliMUONPad::AreOverlapping(pad,tmp,fgkDecreaseSize);
296 //_____________________________________________________________________________
298 AliMUONClusterFinderMLEM::CheckPrecluster(const AliMUONCluster& origCluster)
300 /// Check precluster in order to attempt to simplify it (mostly for
301 /// two-cathode preclusters)
303 if (origCluster.Multiplicity()==1)
305 // Disregard one-pad clusters (leftovers from splitting)
309 Timer(kCheckPreCluster)->Start(kFALSE);
312 AliMUONCluster* cluster = static_cast<AliMUONCluster*>(origCluster.Clone());
316 AliDebug(2,"Start of CheckPreCluster=");
317 StdoutToAliDebug(2,cluster->Print("full"));
319 // Check if one-cathode precluster
320 Int_t i1 = cluster->Multiplicity(0) ? 0 : 1;
321 Int_t i2 = cluster->Multiplicity(1) ? 1 : 0;
323 AliMUONCluster* rv(0x0);
327 rv = CheckPreclusterTwoCathodes(cluster);
331 rv = CheckPreclusterOneCathode(cluster);
333 Timer(kCheckPreCluster)->Stop();
337 //_____________________________________________________________________________
339 AliMUONClusterFinderMLEM::CheckPreclusterOneCathode(AliMUONCluster* cluster)
341 /// Check single-cathode precluster
342 AliWarning("Reimplement me!");
343 AliDebug(2,"End of CheckPreClusterOneCathode=");
344 StdoutToAliDebug(2,cluster->Print("full"));
349 //_____________________________________________________________________________
351 AliMUONClusterFinderMLEM::CheckPreclusterTwoCathodes(AliMUONCluster* cluster)
353 // Check two-cathode cluster
355 Int_t i1 = cluster->Multiplicity(0) ? 0 : 1;
356 Int_t i2 = cluster->Multiplicity(1) ? 1 : 0;
358 Int_t npad = cluster->Multiplicity();
359 Int_t* flags = new Int_t[npad];
360 memset(flags,0,npad*sizeof(Int_t));
362 // Check pad overlaps
363 for ( Int_t i=0; i<npad; ++i)
365 AliMUONPad* padi = cluster->Pad(i);
366 if ( padi->Cathode() != i1 ) continue;
367 for (Int_t j=i+1; j<npad; ++j)
369 AliMUONPad* padj = cluster->Pad(j);
370 if ( padj->Cathode() != i2 ) continue;
371 if ( !AliMUONPad::AreOverlapping(*padi,*padj,fgkDecreaseSize) ) continue;
372 flags[i] = flags[j] = 1; // mark overlapped pads
376 // Check if all pads overlap
378 for (Int_t i=0; i<npad; ++i)
380 if (flags[i]) continue;
386 // not all pads overlap.
387 for (Int_t i=0; i<npad; ++i)
389 AliMUONPad* pad = cluster->Pad(i);
390 if (flags[i]) continue;
391 Int_t cath = pad->Cathode();
392 Int_t cath1 = TMath::Even(cath);
393 // Check for edge effect (missing pads on the _other_ cathode)
394 AliMpPad mpPad = fSegmentation[cath1]->PadByPosition(pad->Position(),kFALSE);
395 if (!mpPad.IsValid()) continue;
396 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
397 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),pad->Charge()));
398 cluster->RemovePad(pad);
399 fPreCluster->Pad(i)->Release();
404 // Check correlations of cathode charges
405 if ( !cluster->IsSaturated() && cluster->ChargeAsymmetry()*2 > 1 )
408 Int_t cathode = cluster->MaxRawChargeCathode();
414 // get min and max pad charges on the cathode opposite to the
415 // max pad (given by MaxRawChargeCathode())
417 for ( Int_t i = 0; i < cluster->Multiplicity(); ++i )
419 AliMUONPad* pad = cluster->Pad(i);
420 if ( pad->Cathode() != cathode || !pad->IsReal() )
422 // only consider pads in the opposite cathode, and
423 // onyl consider real pads (i.e. exclude the virtual ones)
426 if ( pad->Charge() < cmin )
428 cmin = pad->Charge();
431 if ( pad->Charge() > cmax )
433 cmax = pad->Charge();
437 AliDebug(2,Form("Pad imin,imax %d,%d cmin,cmax %e,%e",
438 imin,imax,cmin,cmax));
440 // arrange pads according to their distance to the max, normalized
442 Double_t* dist = new Double_t[cluster->Multiplicity()];
447 AliMUONPad* padmax = cluster->Pad(imax);
449 for ( Int_t i = 0; i < cluster->Multiplicity(); ++i )
452 if ( i == imax) continue;
453 AliMUONPad* pad = cluster->Pad(i);
454 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
455 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
456 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
457 dist[i] = TMath::Sqrt(dx*dx+dy*dy);
460 dmin = dist[i] + 1E-3; // distance to the pad with minimum charge
466 TMath::Sort(cluster->Multiplicity(),dist,flags,kFALSE); // in ascending order
468 TObjArray toBeRemoved;
470 for ( Int_t i = 0; i < cluster->Multiplicity(); ++i )
472 Int_t indx = flags[i];
473 AliMUONPad* pad = cluster->Pad(indx);
474 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
475 if ( dist[indx] > dmin )
477 // farther than the minimum pad
478 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
479 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
482 if (dx >= 0 && dy >= 0) continue;
483 if (TMath::Abs(dx) > TMath::Abs(dy) && dx >= 0) continue;
484 if (TMath::Abs(dy) > TMath::Abs(dx) && dy >= 0) continue;
486 if ( pad->Charge() <= cmax || TMath::Abs(dist[indx]-xmax) < 1E-3 )
489 if (TMath::Abs(dist[indx]-xmax) < 1.e-3)
491 cmax = TMath::Max(pad->Charge(),cmax);
495 cmax = pad->Charge();
498 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
499 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),
502 toBeRemoved.AddLast(pad);
503 fPreCluster->Pad(indx)->Release();
506 for ( Int_t i = 0; i <= toBeRemoved.GetLast(); ++i )
508 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.At(i)));
515 AliDebug(2,"End of CheckPreClusterTwoCathodes=");
516 StdoutToAliDebug(2,cluster->Print("full"));
521 //_____________________________________________________________________________
523 AliMUONClusterFinderMLEM::CheckOverlaps()
525 /// For debug only : check if some pixels overlap...
527 Int_t nPix = fPixArray->GetLast()+1;
530 for ( Int_t i = 0; i < nPix; ++i )
532 AliMUONPad* pixelI = Pixel(i);
533 AliMUONPad pi(dummy,dummy,dummy,dummy,
534 pixelI->Coord(0),pixelI->Coord(1),
535 pixelI->Size(0),pixelI->Size(1),0.0);
537 for ( Int_t j = i+1; j < nPix; ++j )
539 AliMUONPad* pixelJ = Pixel(j);
540 AliMUONPad pj(dummy,dummy,dummy,dummy,
541 pixelJ->Coord(0),pixelJ->Coord(1),
542 pixelJ->Size(0),pixelJ->Size(1),0.0);
545 if ( AliMUONPad::AreOverlapping(pi,pj,fgkDecreaseSize,area) )
547 AliInfo(Form("The following 2 pixels (%d and %d) overlap !",i,j));
548 StdoutToAliInfo(pixelI->Print();
549 cout << " Surface = " << pixelI->Size(0)*pixelI->Size(1)*4 << endl;
551 cout << " Surface = " << pixelJ->Size(0)*pixelJ->Size(1)*4 << endl;
552 cout << " Area surface = " << area.Dimensions().X()*area.Dimensions().Y()*4 << endl;
553 cout << "-------" << endl;
561 //_____________________________________________________________________________
562 void AliMUONClusterFinderMLEM::BuildPixArray(AliMUONCluster& cluster)
564 /// Build pixel array for MLEM method
566 Int_t npad = cluster.Multiplicity();
569 AliWarning("Got no pad at all ?!");
574 if ( !cluster.Multiplicity(0) || !cluster.Multiplicity(1) )
576 BuildPixArrayOneCathode(cluster);
580 BuildPixArrayTwoCathodes(cluster);
583 fPixArray->Sort(); // FIXME : not really needed, only to compare with ClusterFinderAZ
585 Int_t nPix = fPixArray->GetLast()+1;
587 // AliDebug(2,Form("nPix after BuildPixArray=%d",nPix));
589 Double_t xPadMin(1E9);
590 Double_t yPadMin(1E9);
592 for ( Int_t i = 0; i < cluster.Multiplicity(); ++i )
594 AliMUONPad* pad = cluster.Pad(i);
595 xPadMin = TMath::Min (xPadMin, pad->DX());
596 yPadMin = TMath::Min (yPadMin, pad->DY());
602 for ( Int_t i = 0; i < nPix; ++i )
604 AliMUONPad* pixPtr = Pixel(i);
605 wxmin = TMath::Min(wxmin, pixPtr->Size(0));
606 wymin = TMath::Min(wymin, pixPtr->Size(1));
609 wxmin = TMath::Abs (wxmin - xPadMin/2) > 0.001 ? xPadMin : xPadMin / 2;
610 wymin = TMath::Abs (wymin - yPadMin/2) > 0.001 ? yPadMin : yPadMin / 2;
612 // Check if small pixel X-size
613 AdjustPixel(cluster,wxmin, 0);
614 // Check if small pixel Y-size
615 AdjustPixel(cluster,wymin, 1);
616 // Check if large pixel size
617 AdjustPixel(wxmin, wymin);
619 // Remove discarded pixels
620 for (Int_t i=0; i<nPix; ++i)
622 AliMUONPad* pixPtr = Pixel(i);
623 if (pixPtr->Charge() < 1)
625 AliDebug(2,Form("Removing pixel %d with charge<1 : ",i));
626 StdoutToAliDebug(2,pixPtr->Print());
631 fPixArray->Compress();
632 nPix = fPixArray->GetEntriesFast();
634 // AliDebug(2,Form("nPix after AdjustPixel=%d",nPix));
636 if ( nPix > cluster.Multiplicity() )
638 // AliDebug(2,Form("Will trim number of pixels to number of pads"));
640 // Too many pixels - sort and remove pixels with the lowest signal
642 for ( Int_t i = cluster.Multiplicity(); i<nPix; ++i )
646 fPixArray->Compress();
647 nPix = fPixArray->GetEntriesFast();
648 } // if (nPix > npad)
650 // StdoutToAliDebug(2,cout << "End of BuildPixelArray:" << endl;
651 // fPixArray->Print(););
652 CheckOverlaps();//FIXME : this is for debug only. Remove it.
655 //_____________________________________________________________________________
656 void AliMUONClusterFinderMLEM::BuildPixArrayOneCathode(AliMUONCluster& cluster)
658 /// From a single-cathode cluster, build the pixel array
660 // AliDebug(2,Form("cluster.Multiplicity=%d",cluster.Multiplicity()));
662 for ( Int_t j=0; j<cluster.Multiplicity(); ++j)
664 AliMUONPad* pad = cluster.Pad(j);
665 AliMUONPad* pixPtr = new AliMUONPad(pad->Position(),pad->Dimensions(),
667 fPixArray->Add(pixPtr);
671 //_____________________________________________________________________________
672 void AliMUONClusterFinderMLEM::BuildPixArrayTwoCathodes(AliMUONCluster& cluster)
674 /// From a two-cathodes cluster, build the pixel array
676 // AliDebug(2,Form("cluster.Multiplicity=%d",cluster.Multiplicity()));
678 Int_t i1 = cluster.Pad(0)->Cathode();
679 Int_t i2 = TMath::Even(i1);
681 for ( Int_t i = 0; i < cluster.Multiplicity(); ++i)
683 AliMUONPad* padi = cluster.Pad(i);
684 if (padi->Cathode() != i1) continue;
686 for ( Int_t j = 1; j < cluster.Multiplicity(); ++j)
688 AliMUONPad* padj = cluster.Pad(j);
689 if (padj->Cathode() != i2) continue;
693 if ( AliMUONPad::AreOverlapping(*padi,*padj,fgkDecreaseSize,overlap) )
695 AliMUONPad* pixPtr = new AliMUONPad(overlap.Position(),
696 overlap.Dimensions(),
697 TMath::Min(padi->Charge(),padj->Charge()));
698 if ( ( padi->Charge() <= padj->Charge() && padi->IsSaturated() ) ||
699 ( padj->Charge() < padi->Charge() && padj->IsSaturated() ) )
701 // if smallest charge of the 2 pads is already above saturation, then
702 // the pixel is saturated...
703 pixPtr->SetSaturated(kTRUE);
705 pixPtr->SetReal(kFALSE);
706 fPixArray->Add(pixPtr);
712 //_____________________________________________________________________________
713 void AliMUONClusterFinderMLEM::AdjustPixel(AliMUONCluster& cluster,
714 Float_t width, Int_t ixy)
716 /// Check if some pixels have small size (adjust if necessary)
718 AliDebug(2,Form("width=%e ixy=%d",width,ixy));
720 AliMUONPad *pixPtr, *pixPtr1 = 0;
721 Int_t ixy1 = TMath::Even(ixy);
722 Int_t nPix = fPixArray->GetEntriesFast();
724 for (Int_t i=0; i<nPix; i++)
727 if (pixPtr->Charge() < 1) continue; // discarded pixel
728 if (pixPtr->Size(ixy)-width < -1.e-4)
731 for (Int_t j=i+1; j<nPix; j++)
734 if (pixPtr1->Charge() < 1) continue; // discarded pixel
735 if (TMath::Abs(pixPtr1->Size(ixy)-width) < fgkDistancePrecision) continue; // right size
736 if (TMath::Abs(pixPtr1->Coord(ixy1)-pixPtr->Coord(ixy1)) > fgkDistancePrecision) continue; // different rows/columns
737 if (TMath::Abs(pixPtr1->Coord(ixy)-pixPtr->Coord(ixy)) < 2*width)
740 Double_t tmp = pixPtr->Coord(ixy) + pixPtr1->Size(ixy)*
741 TMath::Sign (1., pixPtr1->Coord(ixy) - pixPtr->Coord(ixy));
742 pixPtr->SetCoord(ixy, tmp);
743 pixPtr->SetSize(ixy, width);
744 pixPtr->SetCharge(TMath::Min (pixPtr->Charge(),pixPtr1->Charge()));
745 pixPtr1->SetCharge(0);
749 } // for (Int_t j=i+1;
750 if (pixPtr1 || i == nPix-1)
752 // edge pixel - just increase its size
753 for (Int_t j=0; j<cluster.Multiplicity(); ++j)
755 AliMUONPad* pad = cluster.Pad(j);
756 Double_t d = ( ixy == 0 ) ? pad->X() : ( ixy == 1 ) ? pad->Y() : -1E9;
758 if (pixPtr->Coord(ixy) < d)
760 pixPtr->Shift(ixy, pixPtr->Size(ixy)-width);
764 pixPtr->Shift(ixy, -pixPtr->Size(ixy)+width);
766 pixPtr->SetSize(ixy, width);
770 } // if (pixPtr->Size(ixy)-width < -1.e-4)
771 } // for (Int_t i=0; i<nPix;
774 //_____________________________________________________________________________
775 void AliMUONClusterFinderMLEM::AdjustPixel(Double_t wxmin, Double_t wymin)
777 /// Check if some pixels have large size (adjust if necessary)
779 AliDebug(2,Form("wxmin=%e wymin=%e",wxmin,wymin));
781 Int_t n1[2], n2[2], iOK = 1, nPix = fPixArray->GetEntriesFast();
782 AliMUONPad *pixPtr, pix;
783 Double_t xy0[2] = {9999, 9999}, wxy[2], dist[2] = {0};
785 // Check if large pixel size
786 for (Int_t i = 0; i < nPix; i++) {
787 pixPtr = (AliMUONPad*) fPixArray->UncheckedAt(i);
788 if (pixPtr->Charge() < 1) continue; // discarded pixel
789 if (pixPtr->Size(0) - wxmin < 1.e-4) {
790 if (xy0[0] > 9998) xy0[0] = pixPtr->Coord(0); // position of a "normal" pixel
791 if (pixPtr->Size(1) - wymin < 1.e-4) {
792 if (xy0[1] > 9998) xy0[1] = pixPtr->Coord(1); // position of a "normal" pixel
794 } else iOK = 0; // large pixel
796 iOK = 0; // large pixel
797 if (xy0[1] > 9998 && pixPtr->Size(1) - wymin < 1.e-4) xy0[1] = pixPtr->Coord(1); // "normal" pixel
799 if (xy0[0] < 9998 && xy0[1] < 9998) break;
805 //cout << xy0[0] << " " << xy0[1] << endl;
806 for (Int_t i = 0; i < nPix; i++) {
807 pixPtr = (AliMUONPad*) fPixArray->UncheckedAt(i);
808 if (pixPtr->Charge() < 1) continue; // discarded pixel
811 for (Int_t j = 0; j < 2; j++) {
812 if (pixPtr->Size(j) - wxy[j] < 1.e-4) continue;
813 dist[j] = (pixPtr->Coord(j) - xy0[j]) / wxy[j] / 2; // normalized distance to "normal" pixel
814 n2[j] = TMath::Nint (pixPtr->Size(j) / wxy[j]);
815 n1[j] = n2[j] == 1 ? TMath::Nint(dist[j]) : (Int_t)dist[j];
817 if (n1[0] > 998 && n1[1] > 998) continue;
818 if (fDebug) cout << " Different " << pixPtr->Size(0) << " " << wxy[0] << " "
819 << pixPtr->Size(1) << " " << wxy[1] <<endl;
821 if (n2[0] > 2 || n2[1] > 2) {
822 //cout << n2[0] << " " << n2[1] << endl;
823 if (n2[0] > 2 && n1[0] < 999) n1[0]--;
824 if (n2[1] > 2 && n1[1] < 999) n1[1]--;
826 //cout << n1[0] << " " << n2[0] << " " << n1[1] << " " << n2[1] << endl;
828 pix.SetSize(0, wxy[0]); pix.SetSize(1, wxy[1]);
830 for (Int_t ii = 0; ii < n2[0]; ii++) {
831 if (n1[0] < 999) pix.SetCoord(0, xy0[0] + (n1[0] + TMath::Sign(1.,dist[0]) * ii) * 2 * wxy[0]);
832 for (Int_t jj = 0; jj < n2[1]; jj++) {
833 if (n1[1] < 999) pix.SetCoord(1, xy0[1] + (n1[1] + TMath::Sign(1.,dist[1]) * jj) * 2 * wxy[1]);
834 fPixArray->Add(new AliMUONPad(pix));
838 pixPtr->SetCharge(0);
839 } // for (Int_t i = 0; i < nPix;
842 //_____________________________________________________________________________
844 AliMUONClusterFinderMLEM::Plot(const char* basename)
846 /// Make a plot and save it as png
850 TCanvas* c = new TCanvas("MLEM","MLEM",800,600);
855 c->Print(Form("%s.EVT%d.DE%d.CLU%d.png",basename,fEventNumber,
856 fDetElemId,fClusterNumber));
859 //_____________________________________________________________________________
861 AliMUONClusterFinderMLEM::ComputeCoefficients(AliMUONCluster& cluster,
865 /// Compute coefficients needed for MLEM algorithm
867 Int_t nPix = fPixArray->GetLast()+1;
869 memset(probi,0,nPix*sizeof(Double_t));
871 for ( Int_t j=0; j<cluster.Multiplicity(); ++j )
873 AliMUONPad* pad = cluster.Pad(j);
876 for ( Int_t ipix=0; ipix<nPix; ++ipix )
878 Int_t indx1 = indx + ipix;
879 if (pad->Status() < 0)
884 AliMUONPad* pixPtr = Pixel(ipix);
885 // coef is the charge (given by Mathieson integral) on pad, assuming
886 // the Mathieson is center at pixel.
887 coef[indx1] = fSplitter->ChargeIntegration(pixPtr->Coord(0), pixPtr->Coord(1), *pad);
888 // AliDebug(2,Form("pad=(%d,%d,%e,%e,%e,%e) pix=(%e,%e,%e,%e) coef %e",
889 // pad->Ix(),pad->Iy(),
890 // pad->X(),pad->Y(),
891 // pad->DX(),pad->DY(),
892 // pixPtr->Coord(0),pixPtr->Coord(1),
893 // pixPtr->Size(0),pixPtr->Size(1),
896 probi[ipix] += coef[indx1];
901 //_____________________________________________________________________________
902 Bool_t AliMUONClusterFinderMLEM::MainLoop(AliMUONCluster& cluster, Int_t iSimple)
904 /// Repeat MLEM algorithm until pixel size becomes sufficiently small
906 Int_t nPix = fPixArray->GetLast()+1;
908 AliDebug(2,Form("nPix=%d iSimple=%d, precluster=",nPix,iSimple));
909 StdoutToAliDebug(2,cluster.Print("full"););
913 AliDebug(1,"No pixels, why am I here then ?");
916 AddVirtualPad(cluster); // add virtual pads if necessary
918 Int_t npadTot = cluster.Multiplicity();
920 for (Int_t i = 0; i < npadTot; ++i)
922 if (cluster.Pad(i)->Status() == 0) ++npadOK;
927 Double_t* probi(0x0);
928 Int_t lc(0); // loop counter (for debug)
935 mlem = (TH2D*) gROOT->FindObject("mlem");
938 AliDebug(2,Form("lc %d nPix %d(%d) npadTot %d npadOK %d",lc,nPix,fPixArray->GetLast()+1,npadTot,npadOK));
939 AliDebug(2,Form("EVT%d PixArray=",fEventNumber));
940 StdoutToAliDebug(2,fPixArray->Print("","full"));
942 coef = new Double_t [npadTot*nPix];
943 probi = new Double_t [nPix];
945 // Calculate coefficients and pixel visibilities
946 ComputeCoefficients(cluster,coef,probi);
948 for (Int_t ipix=0; ipix<nPix; ++ipix)
950 if (probi[ipix] < 0.01)
952 AliMUONPad* pixel = Pixel(ipix);
953 AliDebug(2,Form("Setting the following pixel to invisible as its probi<0.01:"));
954 StdoutToAliDebug(2,cout << Form(" -- ipix %3d --- "); pixel->Print(););
955 pixel->SetCharge(0); // "invisible" pixel
960 Mlem(cluster,coef, probi, 15);
962 Double_t xylim[4] = {999, 999, 999, 999};
963 AliMUONPad* pixPtr(0x0);
965 for ( Int_t ipix=0; ipix<nPix; ++ipix )
967 pixPtr = Pixel(ipix);
968 for ( Int_t i=0; i<4; ++i )
970 xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
973 for (Int_t i=0; i<4; i++)
975 xylim[i] -= pixPtr->Size(i/2);
979 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
980 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
982 // StdoutToAliDebug(2,cout << "pixel used for nx,ny computation : "; pixPtr->Print(););
983 AliDebug(2,Form("lc %d pixPtr size = %e,%e nx,ny=%d,%d xylim=%e,%e,%e,%e",
984 lc,pixPtr->Size(0),pixPtr->Size(1),nx,ny,
985 xylim[0],-xylim[1],xylim[2],-xylim[3]
988 mlem = new TH2D("mlem","mlem",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
990 for (Int_t ipix=0; ipix<nPix; ++ipix)
992 AliMUONPad* pixPtr = Pixel(ipix);
993 mlem->Fill(pixPtr->Coord(0),pixPtr->Coord(1),pixPtr->Charge());
996 // Check if the total charge of pixels is too low
998 for ( Int_t i=0; i<nPix; ++i)
1000 qTot += Pixel(i)->Charge();
1003 if ( qTot < 1.e-4 || ( npadOK < 3 && qTot < 7 ) )
1005 AliDebug(1,Form("Deleting the above cluster (charge %e too low, npadOK=%d)",qTot,npadOK));
1010 fPixArray->Delete();
1011 for ( Int_t i=0; i<npadTot; ++i)
1013 AliMUONPad* pad = cluster.Pad(i);
1014 if ( pad->Status() == 0) pad->SetStatus(-1);
1021 // Simple cluster - skip further passes thru EM-procedure
1027 fPixArray->Delete();
1031 // Calculate position of the center-of-gravity around the maximum pixel
1033 FindCOG(mlem, xyCOG);
1035 if (TMath::Min(pixPtr->Size(0),pixPtr->Size(1)) < 0.07 &&
1036 pixPtr->Size(0) > pixPtr->Size(1)) break;
1038 // Sort pixels according to the charge
1040 Double_t pixMin = 0.01*Pixel(0)->Charge();
1041 pixMin = TMath::Min(pixMin,50.);
1043 // Decrease pixel size and shift pixels to make them centered at
1045 Int_t indx = (pixPtr->Size(0)>pixPtr->Size(1)) ? 0 : 1;
1048 Double_t shift[2] = { 0.0, 0.0 };
1049 for (Int_t i=0; i<4; i++) xylim[i] = 999;
1052 for (Int_t ipix=0; ipix<nPix1; ++ipix)
1054 AliMUONPad* pixPtr = Pixel(ipix);
1055 if ( nPix >= npadOK // too many pixels already
1057 pixPtr->Charge() < pixMin // too low charge
1063 for (Int_t i=0; i<2; ++i)
1067 pixPtr->SetCharge(10);
1068 pixPtr->SetSize(indx, pixPtr->Size(indx)/2);
1069 width = -pixPtr->Size(indx);
1070 pixPtr->Shift(indx, width);
1071 // Shift pixel position
1075 for (Int_t j=0; j<2; ++j)
1077 shift[j] = pixPtr->Coord(j) - xyCOG[j];
1078 shift[j] -= ((Int_t)(shift[j]/pixPtr->Size(j)/2))*pixPtr->Size(j)*2;
1081 pixPtr->Shift(0, -shift[0]);
1082 pixPtr->Shift(1, -shift[1]);
1086 pixPtr = new AliMUONPad(*pixPtr);
1087 pixPtr->Shift(indx, -2*width);
1088 fPixArray->Add(pixPtr);
1090 for (Int_t i=0; i<4; ++i)
1092 xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
1094 } // for (Int_t i=0; i<2;
1096 } // for (Int_t ipix=0;
1098 fPixArray->Compress();
1099 nPix = fPixArray->GetEntriesFast();
1101 AliDebug(2,Form("After shift:"));
1102 StdoutToAliDebug(2,fPixArray->Print("","full"););
1103 Plot(Form("mlem.lc%d",lc+1));
1105 AliDebug(2,Form(" xyCOG=%9.6f %9.6f xylim=%9.6f,%9.6f,%9.6f,%9.6f",
1108 xylim[2],xylim[3]));
1110 // Remove excessive pixels
1113 for (Int_t ipix=npadOK; ipix<nPix; ++ipix)
1120 AliMUONPad* pixPtr = Pixel(0);
1121 // add pixels if the maximum is at the limit of pixel area
1122 // start from Y-direction
1124 for (Int_t i=3; i>-1; --i)
1126 if (nPix < npadOK &&
1127 TMath::Abs((i%2 ? -1 : 1)*xylim[i]-xyCOG[i/2]) < pixPtr->Size(i/2))
1129 AliMUONPad* p = static_cast<AliMUONPad*>(pixPtr->Clone());
1130 p->SetCoord(i/2, xyCOG[i/2]+(i%2 ? 2:-2)*pixPtr->Size(i/2));
1131 j = TMath::Even (i/2);
1132 p->SetCoord(j, xyCOG[j]);
1133 AliDebug(2,Form("Adding pixel on the edge (i=%d) ",i));
1134 StdoutToAliDebug(2,cout << " ---- ";
1135 p->Print("corners"););
1141 fPixArray->Compress();
1142 nPix = fPixArray->GetEntriesFast();
1149 AliDebug(2,Form("At the end of while loop nPix=%d : ",fPixArray->GetLast()+1));
1150 StdoutToAliDebug(2,fPixArray->Print("","full"););
1152 // remove pixels with low signal or low visibility
1153 // Cuts are empirical !!!
1154 Double_t thresh = TMath::Max (mlem->GetMaximum()/100.,1.);
1155 thresh = TMath::Min (thresh,50.);
1157 Double_t charge = 0;
1159 for ( Int_t i=0; i<nPix; ++i)
1161 cmax = TMath::Max (cmax,probi[i]);
1164 // Mark pixels which should be removed
1165 for (Int_t i=0; i<nPix; ++i)
1167 AliMUONPad* pixPtr = Pixel(i);
1168 charge = pixPtr->Charge();
1169 if (charge < thresh)
1171 pixPtr->SetCharge(-charge);
1175 // Move charge of removed pixels to their nearest neighbour (to keep total charge the same)
1177 for (Int_t i=0; i<nPix; ++i)
1179 AliMUONPad* pixPtr = Pixel(i);
1180 charge = pixPtr->Charge();
1181 if (charge > 0) continue;
1182 near = FindNearest(pixPtr);
1183 pixPtr->SetCharge(0);
1184 probi[i] = 0; // make it "invisible"
1185 AliMUONPad* pnear = Pixel(near);
1186 pnear->SetCharge(pnear->Charge() + (-charge));
1188 Mlem(cluster,coef,probi,2);
1190 AliDebug(2,Form("Before splitting nPix=%d EVT %d DE %d",fPixArray->GetLast()+1,fEventNumber,fDetElemId));
1191 StdoutToAliDebug(2,fPixArray->Print("","full"););
1192 Plot("mlem.beforesplit");
1195 for (Int_t i=0; i<nPix; ++i)
1197 AliMUONPad* pixPtr = Pixel(i);
1198 Int_t ix = mlem->GetXaxis()->FindBin(pixPtr->Coord(0));
1199 Int_t iy = mlem->GetYaxis()->FindBin(pixPtr->Coord(1));
1200 mlem->SetBinContent(ix, iy, pixPtr->Charge());
1203 // Try to split into clusters
1205 if (mlem->GetSum() < 1)
1211 fSplitter->Split(cluster,mlem,coef,fClusterList);
1218 fPixArray->Delete();
1223 //_____________________________________________________________________________
1224 void AliMUONClusterFinderMLEM::Mlem(AliMUONCluster& cluster,
1225 Double_t* coef, Double_t* probi,
1228 /// Use MLEM to find pixel charges
1230 Int_t nPix = fPixArray->GetEntriesFast();
1232 Int_t npad = cluster.Multiplicity();
1234 Double_t* probi1 = new Double_t[nPix];
1235 Double_t probMax = 0;
1236 Double_t tmp = TMath::MaxElement(nPix,probi);
1238 for (Int_t ipix=0; ipix<nPix; ++ipix)
1240 probMax = TMath::Max(probMax,probi[ipix]);
1243 if (probMax!=tmp) { AliWarning(Form("probMax=%e tmp=%e",probMax,tmp)); }
1245 for (Int_t iter=0; iter<nIter; ++iter)
1248 for (Int_t ipix=0; ipix<nPix; ++ipix)
1250 // Correct each pixel
1251 if (probi[ipix] < 0.01) continue; // skip "invisible" pixel
1253 probi1[ipix] = probMax;
1254 for (Int_t j=0; j<npad; j++)
1256 AliMUONPad* pad = cluster.Pad(j);
1257 if (pad->Status() < 0) continue;
1259 Int_t indx1 = j*nPix;
1260 Int_t indx = indx1 + ipix;
1261 // Calculate expectation
1262 for (Int_t i=0; i<nPix; i++)
1264 sum1 += Pixel(i)->Charge()*coef[indx1+i];
1266 if ( pad->Charge() > fgkSaturation-1 && sum1 > pad->Charge() ) //FIXME : remove usage of fgkSaturation
1268 if ( !pad->IsSaturated() )
1270 AliWarning("Got a pad charge above saturation not backed-up by pad->IsSaturated() function : ");
1271 StdoutToAliWarning(pad->Print("full"));
1273 // correct for pad charge overflows
1274 probi1[ipix] -= coef[indx];
1278 if (coef[indx] > 1.e-6)
1280 sum += pad->Charge()*coef[indx]/sum1;
1282 } // for (Int_t j=0;
1283 AliMUONPad* pixPtr = Pixel(ipix);
1284 if (probi1[ipix] > 1.e-6)
1286 pixPtr->SetCharge(pixPtr->Charge()*sum/probi1[ipix]);
1288 } // for (Int_t ipix=0;
1289 } // for (Int_t iter=0;
1293 //_____________________________________________________________________________
1294 void AliMUONClusterFinderMLEM::FindCOG(TH2D *mlem, Double_t *xyc)
1296 /// Calculate position of the center-of-gravity around the maximum pixel
1298 Int_t ixmax, iymax, ix, nsumx=0, nsumy=0, nsum=0;
1299 Int_t i1 = -9, j1 = -9;
1300 mlem->GetMaximumBin(ixmax,iymax,ix);
1301 Int_t nx = mlem->GetNbinsX();
1302 Int_t ny = mlem->GetNbinsY();
1303 Double_t thresh = mlem->GetMaximum()/10;
1304 Double_t x, y, cont, xq=0, yq=0, qq=0;
1306 for (Int_t i=TMath::Max(1,iymax-1); i<=TMath::Min(ny,iymax+1); i++) {
1307 y = mlem->GetYaxis()->GetBinCenter(i);
1308 for (Int_t j=TMath::Max(1,ixmax-1); j<=TMath::Min(nx,ixmax+1); j++) {
1309 cont = mlem->GetCellContent(j,i);
1310 if (cont < thresh) continue;
1311 if (i != i1) {i1 = i; nsumy++;}
1312 if (j != j1) {j1 = j; nsumx++;}
1313 x = mlem->GetXaxis()->GetBinCenter(j);
1322 Int_t i2 = 0, j2 = 0;
1325 // one bin in Y - add one more (with the largest signal)
1326 for (Int_t i=TMath::Max(1,iymax-1); i<=TMath::Min(ny,iymax+1); i++) {
1327 if (i == iymax) continue;
1328 for (Int_t j=TMath::Max(1,ixmax-1); j<=TMath::Min(nx,ixmax+1); j++) {
1329 cont = mlem->GetCellContent(j,i);
1332 x = mlem->GetXaxis()->GetBinCenter(j);
1333 y = mlem->GetYaxis()->GetBinCenter(i);
1342 if (i2 != i1) nsumy++;
1343 if (j2 != j1) nsumx++;
1345 } // if (nsumy == 1)
1348 // one bin in X - add one more (with the largest signal)
1350 for (Int_t j=TMath::Max(1,ixmax-1); j<=TMath::Min(nx,ixmax+1); j++) {
1351 if (j == ixmax) continue;
1352 for (Int_t i=TMath::Max(1,iymax-1); i<=TMath::Min(ny,iymax+1); i++) {
1353 cont = mlem->GetCellContent(j,i);
1356 x = mlem->GetXaxis()->GetBinCenter(j);
1357 y = mlem->GetYaxis()->GetBinCenter(i);
1366 if (i2 != i1) nsumy++;
1367 if (j2 != j1) nsumx++;
1369 } // if (nsumx == 1)
1371 xyc[0] = xq/qq; xyc[1] = yq/qq;
1372 AliDebug(2,Form("x,y COG = %e,%e",xyc[0],xyc[1]));
1375 //_____________________________________________________________________________
1376 Int_t AliMUONClusterFinderMLEM::FindNearest(AliMUONPad *pixPtr0)
1378 /// Find the pixel nearest to the given one
1379 /// (algorithm may be not very efficient)
1381 Int_t nPix = fPixArray->GetEntriesFast(), imin = 0;
1382 Double_t rmin = 99999, dx = 0, dy = 0, r = 0;
1383 Double_t xc = pixPtr0->Coord(0), yc = pixPtr0->Coord(1);
1386 for (Int_t i=0; i<nPix; i++) {
1387 pixPtr = (AliMUONPad*) fPixArray->UncheckedAt(i);
1388 if (pixPtr->Charge() < 0.5) continue;
1389 dx = (xc - pixPtr->Coord(0)) / pixPtr->Size(0);
1390 dy = (yc - pixPtr->Coord(1)) / pixPtr->Size(1);
1391 r = dx *dx + dy * dy;
1392 if (r < rmin) { rmin = r; imin = i; }
1398 //_____________________________________________________________________________
1400 AliMUONClusterFinderMLEM::Timer(Int_t i) const
1402 /// Return timer at index i
1403 return static_cast<TStopwatch*>(fTimers->At(i));
1406 //_____________________________________________________________________________
1408 AliMUONClusterFinderMLEM::Paint(Option_t*)
1410 /// Paint cluster and pixels
1412 AliMpArea area(fPreCluster->Area());
1414 gPad->Range(area.LeftBorder(),area.DownBorder(),area.RightBorder(),area.UpBorder());
1416 gVirtualX->SetFillStyle(1001);
1417 gVirtualX->SetFillColor(3);
1418 gVirtualX->SetLineColor(3);
1422 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1424 AliMUONPad* pixel = Pixel(i);
1426 gPad->PaintBox(pixel->Coord(0)-pixel->Size(0)*s,
1427 pixel->Coord(1)-pixel->Size(1)*s,
1428 pixel->Coord(0)+pixel->Size(0)*s,
1429 pixel->Coord(1)+pixel->Size(1)*s);
1431 // for ( Int_t sign = -1; sign < 2; sign +=2 )
1433 // gPad->PaintLine(pixel->Coord(0) - pixel->Size(0),
1434 // pixel->Coord(1) + sign*pixel->Size(1),
1435 // pixel->Coord(0) + pixel->Size(0),
1436 // pixel->Coord(1) - sign*pixel->Size(1)
1442 gVirtualX->SetFillStyle(0);
1444 fPreCluster->Paint();
1446 gVirtualX->SetLineColor(1);
1447 gVirtualX->SetLineWidth(2);
1448 gVirtualX->SetFillStyle(0);
1449 gVirtualX->SetTextColor(1);
1450 gVirtualX->SetTextAlign(22);
1452 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1454 AliMUONPad* pixel = Pixel(i);
1455 gPad->PaintBox(pixel->Coord(0)-pixel->Size(0),
1456 pixel->Coord(1)-pixel->Size(1),
1457 pixel->Coord(0)+pixel->Size(0),
1458 pixel->Coord(1)+pixel->Size(1));
1459 gVirtualX->SetTextSize(pixel->Size(0)*60);
1461 gPad->PaintText(pixel->Coord(0),pixel->Coord(1),Form("%d",(Int_t)(pixel->Charge())));
1465 //_____________________________________________________________________________
1466 Int_t AliMUONClusterFinderMLEM::FindLocalMaxima(TObjArray *pixArray, Int_t *localMax, Double_t *maxVal)
1468 /// Find local maxima in pixel space for large preclusters in order to
1469 /// try to split them into smaller pieces (to speed up the MLEM procedure)
1470 /// or to find additional fitting seeds if clusters were not completely resolved
1472 AliDebug(1,Form("nPix=%d",pixArray->GetLast()+1));
1475 //if (pixArray == fPixArray) hist = (TH2D*) gROOT->FindObject("anode");
1476 //else { hist = (TH2D*) gROOT->FindObject("anode1"); cout << hist << endl; }
1477 //if (hist) hist->Delete();
1479 Double_t xylim[4] = {999, 999, 999, 999};
1481 Int_t nPix = pixArray->GetEntriesFast();
1482 AliMUONPad *pixPtr = 0;
1483 for (Int_t ipix=0; ipix<nPix; ipix++) {
1484 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
1485 for (Int_t i=0; i<4; i++)
1486 xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
1488 for (Int_t i=0; i<4; i++) xylim[i] -= pixPtr->Size(i/2);
1490 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
1491 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
1492 if (pixArray == fPixArray) hist = new TH2D("anode","anode",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1493 else hist = new TH2D("anode1","anode1",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1494 for (Int_t ipix=0; ipix<nPix; ipix++) {
1495 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
1496 hist->Fill(pixPtr->Coord(0), pixPtr->Coord(1), pixPtr->Charge());
1498 // if (fDraw && pixArray == fPixArray) fDraw->DrawHist("c2", hist);
1500 Int_t nMax = 0, indx;
1501 Int_t *isLocalMax = new Int_t[ny*nx];
1502 for (Int_t i=0; i<ny*nx; i++) isLocalMax[i] = 0;
1504 for (Int_t i=1; i<=ny; i++) {
1506 for (Int_t j=1; j<=nx; j++) {
1507 if (hist->GetCellContent(j,i) < 0.5) continue;
1508 //if (isLocalMax[indx+j-1] < 0) continue;
1509 if (isLocalMax[indx+j-1] != 0) continue;
1510 FlagLocalMax(hist, i, j, isLocalMax);
1514 for (Int_t i=1; i<=ny; i++) {
1516 for (Int_t j=1; j<=nx; j++) {
1517 if (isLocalMax[indx+j-1] > 0) {
1518 localMax[nMax] = indx + j - 1;
1519 maxVal[nMax++] = hist->GetCellContent(j,i);
1520 if (nMax > 99) AliFatal(" Too many local maxima !!!");
1524 if (fDebug) cout << " Local max: " << nMax << endl;
1525 delete [] isLocalMax; isLocalMax = 0;
1529 //_____________________________________________________________________________
1530 void AliMUONClusterFinderMLEM::FlagLocalMax(TH2D *hist, Int_t i, Int_t j, Int_t *isLocalMax)
1532 /// Flag pixels (whether or not local maxima)
1534 Int_t nx = hist->GetNbinsX();
1535 Int_t ny = hist->GetNbinsY();
1536 Int_t cont = TMath::Nint (hist->GetCellContent(j,i));
1537 Int_t cont1 = 0, indx = (i-1)*nx+j-1, indx1 = 0, indx2 = 0;
1539 for (Int_t i1=i-1; i1<i+2; i1++) {
1540 if (i1 < 1 || i1 > ny) continue;
1541 indx1 = (i1 - 1) * nx;
1542 for (Int_t j1=j-1; j1<j+2; j1++) {
1543 if (j1 < 1 || j1 > nx) continue;
1544 if (i == i1 && j == j1) continue;
1545 indx2 = indx1 + j1 - 1;
1546 cont1 = TMath::Nint (hist->GetCellContent(j1,i1));
1547 if (cont < cont1) { isLocalMax[indx] = -1; return; }
1548 else if (cont > cont1) isLocalMax[indx2] = -1;
1549 else { // the same charge
1550 isLocalMax[indx] = 1;
1551 if (isLocalMax[indx2] == 0) {
1552 FlagLocalMax(hist, i1, j1, isLocalMax);
1553 if (isLocalMax[indx2] < 0) { isLocalMax[indx] = -1; return; }
1554 else isLocalMax[indx2] = -1;
1559 isLocalMax[indx] = 1; // local maximum
1562 //_____________________________________________________________________________
1563 void AliMUONClusterFinderMLEM::FindCluster(AliMUONCluster& cluster,
1564 Int_t *localMax, Int_t iMax)
1566 /// Find pixel cluster around local maximum \a iMax and pick up pads
1567 /// overlapping with it
1569 TH2D *hist = (TH2D*) gROOT->FindObject("anode");
1570 Int_t nx = hist->GetNbinsX();
1571 Int_t ny = hist->GetNbinsY();
1572 Int_t ic = localMax[iMax] / nx + 1;
1573 Int_t jc = localMax[iMax] % nx + 1;
1574 Bool_t *used = new Bool_t[ny*nx];
1575 for (Int_t i=0; i<ny*nx; i++) used[i] = kFALSE;
1577 // Drop all pixels from the array - pick up only the ones from the cluster
1578 fPixArray->Delete();
1580 Double_t wx = hist->GetXaxis()->GetBinWidth(1)/2;
1581 Double_t wy = hist->GetYaxis()->GetBinWidth(1)/2;
1582 Double_t yc = hist->GetYaxis()->GetBinCenter(ic);
1583 Double_t xc = hist->GetXaxis()->GetBinCenter(jc);
1584 Double_t cont = hist->GetCellContent(jc,ic);
1585 fPixArray->Add(new AliMUONPad (xc, yc, wx, wy, cont));
1586 used[(ic-1)*nx+jc-1] = kTRUE;
1587 fSplitter->AddBin(hist, ic, jc, 1, used, (TObjArray*)0); // recursive call
1589 Int_t nPix = fPixArray->GetEntriesFast();
1590 Int_t npad = cluster.Multiplicity();
1592 for (Int_t i=0; i<nPix; ++i)
1594 AliMUONPad* pixPtr = Pixel(i);
1595 pixPtr->SetSize(0,wx);
1596 pixPtr->SetSize(1,wy);
1599 // Pick up pads which overlap with found pixels
1600 for (Int_t i=0; i<npad; i++)
1602 cluster.Pad(i)->SetStatus(-1);
1605 for (Int_t i=0; i<nPix; i++)
1607 AliMUONPad* pixPtr = Pixel(i);
1608 for (Int_t j=0; j<npad; ++j)
1610 AliMUONPad* pad = cluster.Pad(j);
1611 if ( Overlap(*pad,*pixPtr) )
1618 delete [] used; used = 0;
1621 //_____________________________________________________________________________
1622 AliMUONClusterFinderMLEM&
1623 AliMUONClusterFinderMLEM::operator=(const AliMUONClusterFinderMLEM& rhs)
1625 /// Protected assignement operator
1627 if (this == &rhs) return *this;
1629 AliFatal("Not implemented.");
1634 //_____________________________________________________________________________
1636 AliMUONClusterFinderMLEM::Neighbours(Int_t cathode, Int_t ix, Int_t iy,
1637 Int_t& n, Int_t* xList, Int_t* yList)
1639 /// Get the list of neighbours of pad at (cathode,ix,iy)
1642 const AliMpVSegmentation* seg = fSegmentation[cathode];
1644 AliMpPad pad = seg->PadByIndices(AliMpIntPair(ix,iy),kTRUE);
1646 // Define the region to look into : a region slightly bigger
1647 // than the pad itself (5% bigger), in order to catch first neighbours.
1649 AliMpArea area(pad.Position(),pad.Dimensions()*1.05);
1651 AliMpVPadIterator* it = seg->CreateIterator(area);
1653 while ( !it->IsDone() && n < 10 )
1655 AliMpPad p = it->CurrentItem();
1656 if ( p != pad ) // skip self
1658 xList[n] = p.GetIndices().GetFirst();
1659 yList[n] = p.GetIndices().GetSecond();
1667 //_____________________________________________________________________________
1668 void AliMUONClusterFinderMLEM::AddVirtualPad(AliMUONCluster& cluster)
1670 /// Add virtual pad (with small charge) to improve fit for some
1671 /// clusters (when pad with max charge is at the extreme of the cluster)
1673 // Get number of pads in X and Y-directions
1674 Int_t nInX = -1, nInY;
1675 PadsInXandY(cluster,nInX, nInY);
1677 if (fDebug) cout << " Chamber: " << fDetElemId / 100 - 1 << " " << nInX << " " << nInY << endl;
1679 // Add virtual pad only if number of pads per direction == 2
1680 if (nInX != 2 && nInY != 2) return;
1684 // Find pads with max charge
1685 Int_t maxpad[2][2] = {{-1, -1}, {-1, -1}}, cath;
1686 Double_t sigmax[2] = {0}, aamax[2] = {0};
1687 for (Int_t j=0; j<cluster.Multiplicity(); ++j)
1689 AliMUONPad* pad = cluster.Pad(j);
1690 if (pad->Status() != 0) continue;
1691 cath = pad->Cathode();
1692 if (pad->Charge() > sigmax[cath])
1694 maxpad[cath][1] = maxpad[cath][0];
1695 aamax[cath] = sigmax[cath];
1696 sigmax[cath] = pad->Charge();
1697 maxpad[cath][0] = j;
1701 if (maxpad[0][0] >= 0 && maxpad[0][1] < 0 || maxpad[1][0] >= 0 && maxpad[1][1] < 0)
1703 for (Int_t j=0; j<cluster.Multiplicity(); ++j)
1705 AliMUONPad* pad = cluster.Pad(j);
1706 if (pad->Status() != 0) continue;
1707 cath = pad->Cathode();
1708 if (j == maxpad[cath][0] || j == maxpad[cath][1]) continue;
1709 if ( pad->Charge() > aamax[cath])
1711 aamax[cath] = pad->Charge();
1712 maxpad[cath][1] = j;
1717 // cout << "-------AddVirtualPad" << endl;
1718 // cout << Form("nInX %2d nInY %2d",nInX,nInY) << endl;
1720 // cluster.Print("full");
1722 // for ( Int_t i = 0; i < 2; ++i )
1724 // for ( Int_t j = 0; j < 2; ++j )
1726 // cout << Form("maxpad[%d][%d]=%d",i,j,maxpad[i][j]) << endl;
1730 // Check for mirrors (side X on cathode 0)
1731 Bool_t mirror = kFALSE;
1732 if (maxpad[0][0] >= 0 && maxpad[1][0] >= 0)
1734 AliMUONPad* maxPadCath[] = { cluster.Pad(maxpad[0][0]), cluster.Pad(maxpad[1][0]) };
1735 mirror = maxPadCath[0]->DX() < maxPadCath[0]->DY();
1736 if (!mirror && TMath::Abs( maxPadCath[0]->DX() - maxPadCath[1]->DX()) < 0.001)
1738 // Special case when pads on both cathodes have the same size
1740 for (Int_t j = 0; j < cluster.Multiplicity(); ++j)
1742 AliMUONPad* pad = cluster.Pad(j);
1743 cath = pad->Cathode();
1744 if (j == maxpad[cath][0]) continue;
1745 if ( pad->Ix() != maxPadCath[cath]->Ix() ) continue;
1746 if ( TMath::Abs(pad->Iy() - maxPadCath[cath]->Iy()) == 1 )
1751 if (!yud[0]) mirror = kTRUE; // take the other cathode
1752 } // if (!mirror &&...
1753 } // if (maxpad[0][0] >= 0 && maxpad[1][0] >= 0)
1755 // // Find neughbours of pads with max charges
1756 Int_t nn, xList[10], yList[10], ix0, iy0, ix, iy, neighb;
1757 for (cath=0; cath<2; cath++)
1759 if (!cath && maxpad[0][0] < 0) continue; // one-sided cluster - cathode 1
1760 if (cath && maxpad[1][0] < 0) break; // one-sided cluster - cathode 0
1761 if (maxpad[1][0] >= 0)
1765 if (!cath && nInY != 2) continue;
1766 if (cath && nInX != 2 && (maxpad[0][0] >= 0 || nInY != 2)) continue;
1770 if (!cath && nInX != 2) continue;
1771 if (cath && nInY != 2 && (maxpad[0][0] >= 0 || nInX != 2)) continue;
1775 Int_t iAddX = 0, iAddY = 0, ix1 = 0, iy1 = 0, iPad = 0;
1776 if (maxpad[0][0] < 0) iPad = 1;
1778 for (iPad=0; iPad<2; iPad++)
1780 if (maxpad[cath][iPad] < 0) continue;
1781 if (iPad && !iAddX && !iAddY) break;
1782 if (iPad && cluster.Pad(maxpad[cath][1])->Charge() / sigmax[cath] < 0.5) break;
1784 Int_t neighbx = 0, neighby = 0;
1785 ix0 = cluster.Pad(maxpad[cath][iPad])->Ix();
1786 iy0 = cluster.Pad(maxpad[cath][iPad])->Iy();
1787 Neighbours(cath,ix0,iy0,nn,xList,yList);
1789 for (Int_t j=0; j<nn; j++) {
1790 if (TMath::Abs(xList[j]-ix0) == 1 || xList[j]*ix0 == -1) neighbx++;
1791 if (TMath::Abs(yList[j]-iy0) == 1 || yList[j]*iy0 == -1) neighby++;
1794 if (cath) neighb = neighbx;
1795 else neighb = neighby;
1796 if (maxpad[0][0] < 0) neighb += neighby;
1797 else if (maxpad[1][0] < 0) neighb += neighbx;
1799 if (!cath) neighb = neighbx;
1800 else neighb = neighby;
1801 if (maxpad[0][0] < 0) neighb += neighbx;
1802 else if (maxpad[1][0] < 0) neighb += neighby;
1805 for (Int_t j=0; j< cluster.Multiplicity(); ++j)
1807 AliMUONPad* pad = cluster.Pad(j);
1808 if ( pad->Cathode() != cath) continue;
1811 if (iy == iy0 && ix == ix0) continue;
1812 for (Int_t k=0; k<nn; ++k)
1814 if (xList[k] != ix || yList[k] != iy) continue;
1817 if ((!cath || maxpad[0][0] < 0) &&
1818 (TMath::Abs(iy-iy0) == 1 || iy*iy0 == -1)) {
1819 if (!iPad && TMath::Abs(ix-ix0) == 1 || ix*ix0 == -1) ix1 = xList[k]; //19-12-05
1820 xList[k] = yList[k] = 0;
1824 if ((cath || maxpad[1][0] < 0) &&
1825 (TMath::Abs(ix-ix0) == 1 || ix*ix0 == -1)) {
1826 if (!iPad) ix1 = xList[k]; //19-12-05
1827 xList[k] = yList[k] = 0;
1831 if ((!cath || maxpad[0][0] < 0) &&
1832 (TMath::Abs(ix-ix0) == 1 || ix*ix0 == -1)) {
1833 if (!iPad) ix1 = xList[k]; //19-12-05
1834 xList[k] = yList[k] = 0;
1838 if ((cath || maxpad[1][0] < 0) &&
1839 (TMath::Abs(iy-iy0) == 1 || iy*iy0 == -1)) {
1840 xList[k] = yList[k] = 0;
1845 } // for (Int_t k=0; k<nn;
1847 } // for (Int_t j=0; j< cluster.Multiplicity();
1848 if (!neighb) continue;
1853 for (Int_t j=0; j<nn; j++)
1855 if (xList[j] == 0 && yList[j] == 0) continue;
1856 // npads = fnPads[0] + fnPads[1];
1857 // fPadIJ[0][npads] = cath;
1858 // fPadIJ[1][npads] = 0;
1861 if (TMath::Abs(ix-ix0) == 1 || ix*ix0 == -1) {
1862 if (iy != iy0) continue; // new segmentation - check
1863 if (nInX != 2) continue; // new
1865 if (!cath && maxpad[1][0] >= 0) continue;
1867 if (cath && maxpad[0][0] >= 0) continue;
1869 if (iPad && !iAddX) continue;
1870 AliMpPad pad = fSegmentation[cath]->PadByIndices(AliMpIntPair(ix,iy),kTRUE);
1871 // fXyq[0][npads] = pad.Position().X();
1872 // fXyq[1][npads] = pad.Position().Y();
1873 AliMUONPad muonPad(fDetElemId, cath, ix, iy, pad.Position().X(), pad.Position().Y(), 0, 0, 0);
1874 // fSegmentation[cath]->GetPadC(ix, iy, fXyq[0][npads], fXyq[1][npads], zpad);
1875 // if (fXyq[0][npads] > 1.e+5) continue; // temporary fix
1876 if (muonPad.Coord(0) > 1.e+5) continue; // temporary fix
1877 if (ix == ix1) continue; //19-12-05
1878 if (ix1 == ix0) continue;
1879 if (maxpad[1][0] < 0 || mirror && maxpad[0][0] >= 0) {
1880 // if (!iPad) fXyq[2][npads] = TMath::Min (sigmax[0]/100, 5.);
1881 // else fXyq[2][npads] = TMath::Min (aamax[0]/100, 5.);
1882 if (!iPad) muonPad.SetCharge(TMath::Min (sigmax[0]/100, 5.));
1883 else muonPad.SetCharge(TMath::Min (aamax[0]/100, 5.));
1886 // if (!iPad) fXyq[2][npads] = TMath::Min (sigmax[1]/100, 5.);
1887 // else fXyq[2][npads] = TMath::Min (aamax[1]/100, 5.);
1888 if (!iPad) muonPad.SetCharge(TMath::Min (sigmax[1]/100, 5.));
1889 else muonPad.SetCharge(TMath::Min (aamax[1]/100, 5.));
1891 // fXyq[2][npads] = TMath::Max (fXyq[2][npads], (float)1);
1892 if (muonPad.Charge() < 1.) muonPad.SetCharge(1.);
1893 // fXyq[3][npads] = -2; // flag
1894 // fPadIJ[2][npads] = ix;
1895 // fPadIJ[3][npads] = iy;
1896 muonPad.SetSize(0,-2.); //flag
1900 //AliDebug(1,Form("Add virtual pad in X %f %f %f %3d %3d \n",
1901 // fXyq[2][npads], fXyq[0][npads], fXyq[1][npads], ix, iy));
1902 //muonPad.Charge(), muonPad.Coord(0), muonPad.Coord(1), ix, iy));
1903 if (fDebug) printf(" ***** Add virtual pad in X ***** %f %f %f %3d %3d \n",
1904 muonPad.Charge(), muonPad.Coord(0), muonPad.Coord(1), ix, iy);
1905 cluster.AddPad(muonPad); // add pad to the cluster
1909 if (nInY != 2) continue;
1910 if (!mirror && cath && maxpad[0][0] >= 0) continue;
1911 if (mirror && !cath && maxpad[1][0] >= 0) continue;
1912 if (TMath::Abs(iy-iy0) == 1 || TMath::Abs(iy*iy0) == 1) {
1913 if (ix != ix0) continue; // new segmentation - check
1914 if (iPad && !iAddY) continue;
1915 AliMpPad pad = fSegmentation[cath]->PadByIndices(AliMpIntPair(ix,iy),kTRUE);
1916 // fXyq[0][npads] = pad.Position().X();
1917 // fXyq[1][npads] = pad.Position().Y();
1918 // fSegmentation[cath]->GetPadC(ix, iy, fXyq[0][npads], fXyq[1][npads], zpad);
1919 AliMUONPad muonPad(fDetElemId, cath, ix, iy, pad.Position().X(), pad.Position().Y(), 0, 0, 0);
1920 if (iy1 == iy0) continue;
1921 //if (iPad && iy1 == iy0) continue;
1922 if (maxpad[0][0] < 0 || mirror && maxpad[1][0] >= 0) {
1923 // if (!iPad) fXyq[2][npads] = TMath::Min (sigmax[1]/15, fgkZeroSuppression);
1924 // else fXyq[2][npads] = TMath::Min (aamax[1]/15, fgkZeroSuppression);
1925 if (!iPad) muonPad.SetCharge(TMath::Min (sigmax[1]/15, fgkZeroSuppression));
1926 else muonPad.SetCharge(TMath::Min (aamax[1]/15, fgkZeroSuppression));
1929 // if (!iPad) fXyq[2][npads] = TMath::Min (sigmax[0]/15, fgkZeroSuppression);
1930 // else fXyq[2][npads] = TMath::Min (aamax[0]/15, fgkZeroSuppression);
1931 if (!iPad) muonPad.SetCharge(TMath::Min (sigmax[0]/15, fgkZeroSuppression));
1932 else muonPad.SetCharge(TMath::Min (aamax[0]/15, fgkZeroSuppression));
1934 // fXyq[2][npads] = TMath::Max (fXyq[2][npads], (float)1);
1935 if (muonPad.Charge() < 1.) muonPad.SetCharge(1.);
1936 // fXyq[3][npads] = -2; // flag
1937 // fPadIJ[2][npads] = ix;
1938 // fPadIJ[3][npads] = iy;
1939 muonPad.SetSize(0,-2.); //flag
1943 if (fDebug) printf(" ***** Add virtual pad in Y ***** %f %f %f %3d %3d \n",
1944 muonPad.Charge(), muonPad.Coord(0), muonPad.Coord(1), ix, iy);
1945 cluster.AddPad(muonPad); // add pad to the cluster
1948 } // for (Int_t j=0; j<nn;
1949 } // for (Int_t iPad=0;
1950 } // for (cath=0; cath<2;
1953 //_____________________________________________________________________________
1954 void AliMUONClusterFinderMLEM::PadsInXandY(AliMUONCluster& cluster,
1955 Int_t &nInX, Int_t &nInY) const
1957 /// Find number of pads in X and Y-directions (excluding virtual ones and
1960 Int_t statusToTest = 1;
1962 if ( nInX < 0 ) statusToTest = 0;
1964 Bool_t mustMatch(kTRUE);
1966 AliMpIntPair cn = cluster.NofPads(statusToTest,mustMatch);
1968 nInX = cn.GetFirst();
1969 nInY = cn.GetSecond();
1972 //_____________________________________________________________________________
1973 void AliMUONClusterFinderMLEM::RemovePixel(Int_t i)
1975 /// Remove pixel at index i
1976 AliMUONPad* pixPtr = Pixel(i);
1977 fPixArray->RemoveAt(i);
1981 //_____________________________________________________________________________
1982 void AliMUONClusterFinderMLEM::Simple(AliMUONCluster& cluster)
1984 /// Process simple cluster (small number of pads) without EM-procedure
1986 Int_t nForFit = 1, clustFit[1] = {0}, nfit;
1987 Double_t parOk[3] = {0.};
1988 TObjArray *clusters[1];
1989 clusters[0] = fPixArray;
1991 AliDebug(1,Form("nPix=%d",fPixArray->GetLast()+1));
1993 for (Int_t i = 0; i < cluster.Multiplicity(); ++i)
1995 AliMUONPad* pad = cluster.Pad(i);
1996 if ( pad->Charge() > fgkSaturation-1) //FIXME : remove usage of fgkSaturation
2005 nfit = fSplitter->Fit(cluster,1, nForFit, clustFit, clusters, parOk, fClusterList);
2008 //_____________________________________________________________________________
2010 AliMUONClusterFinderMLEM::Pixel(Int_t i) const
2012 /// Returns pixel at index i
2013 return static_cast<AliMUONPad*>(fPixArray->UncheckedAt(i));
2016 //_____________________________________________________________________________
2018 AliMUONClusterFinderMLEM::Print(Option_t* what) const
2021 TString swhat(what);
2023 if ( swhat.Contains("precluster") )
2025 if ( fPreCluster) fPreCluster->Print();