]> git.uio.no Git - u/mrichter/AliRoot.git/blame - MUON/AliMUONClusterFinderMLEM.cxx
Changed interface to allow combined tracking (Laurent)
[u/mrichter/AliRoot.git] / MUON / AliMUONClusterFinderMLEM.cxx
CommitLineData
c0a16418 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
16/* $Id$ */
17
3d1463c8 18//-----------------------------------------------------------------------------
c0a16418 19/// \class AliMUONClusterFinderMLEM
20///
21/// Clusterizer class based on the Expectation-Maximization algorithm
22///
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).
26///
27/// \author Laurent Aphecetche (for the "new" C++ structure) and
28/// Alexander Zinchenko, JINR Dubna, for the hardcore of it ;-)
3d1463c8 29//-----------------------------------------------------------------------------
c0a16418 30
9edefa04 31#include "AliMUONClusterFinderMLEM.h"
c0a16418 32#include "AliLog.h"
33#include "AliMUONCluster.h"
34#include "AliMUONClusterSplitterMLEM.h"
3887d11e 35#include "AliMUONVDigit.h"
c0a16418 36#include "AliMUONPad.h"
37#include "AliMUONPreClusterFinder.h"
38#include "AliMpPad.h"
39#include "AliMpVPadIterator.h"
40#include "AliMpVSegmentation.h"
41#include "AliRunLoader.h"
3887d11e 42#include "AliMUONVDigitStore.h"
a6435fb4 43#include <Riostream.h>
44#include <TH2.h>
45#include <TMinuit.h>
46#include <TCanvas.h>
a6435fb4 47#include <TMath.h>
a1e17193 48#include <TROOT.h>
f6c291ef 49#include "AliCodeTimer.h"
a6435fb4 50
c0a16418 51/// \cond CLASSIMP
52ClassImp(AliMUONClusterFinderMLEM)
53/// \endcond
54
55const Double_t AliMUONClusterFinderMLEM::fgkZeroSuppression = 6; // average zero suppression value
2abdae6e 56//const Double_t AliMUONClusterFinderMLEM::fgkDistancePrecision = 1e-6; // (cm) used to check overlaps and so on
57const Double_t AliMUONClusterFinderMLEM::fgkDistancePrecision = 1e-3; // (cm) used to check overlaps and so on
c0a16418 58const TVector2 AliMUONClusterFinderMLEM::fgkIncreaseSize(-AliMUONClusterFinderMLEM::fgkDistancePrecision,-AliMUONClusterFinderMLEM::fgkDistancePrecision);
59const TVector2 AliMUONClusterFinderMLEM::fgkDecreaseSize(AliMUONClusterFinderMLEM::fgkDistancePrecision,AliMUONClusterFinderMLEM::fgkDistancePrecision);
60
61 TMinuit* AliMUONClusterFinderMLEM::fgMinuit = 0x0;
8c718ff8 62// Status flags for pads
63const Int_t AliMUONClusterFinderMLEM::fgkZero = 0x0; ///< pad "basic" state
64const Int_t AliMUONClusterFinderMLEM::fgkMustKeep = 0x1; ///< do not kill (for pixels)
65const Int_t AliMUONClusterFinderMLEM::fgkUseForFit = 0x10; ///< should be used for fit
66const Int_t AliMUONClusterFinderMLEM::fgkOver = 0x100; ///< processing is over
67const Int_t AliMUONClusterFinderMLEM::fgkModified = 0x1000; ///< modified pad charge
68const Int_t AliMUONClusterFinderMLEM::fgkCoupled = 0x10000; ///< coupled pad
c0a16418 69
70//_____________________________________________________________________________
b1a19e07 71AliMUONClusterFinderMLEM::AliMUONClusterFinderMLEM(Bool_t plot, AliMUONVClusterFinder* clusterFinder)
c0a16418 72 : AliMUONVClusterFinder(),
b1a19e07 73fPreClusterFinder(clusterFinder),
c0a16418 74fPreCluster(0x0),
75fClusterList(),
76fEventNumber(0),
77fDetElemId(-1),
78fClusterNumber(0),
c0a16418 79fCathBeg(0),
80fPixArray(new TObjArray(20)),
9bbd7f60 81fDebug(0),
c0a16418 82fPlot(plot),
c0a16418 83fSplitter(0x0),
84fNClusters(0),
85fNAddVirtualPads(0)
86{
87 /// Constructor
88
89 fSegmentation[1] = fSegmentation[0] = 0x0;
90
91 fPadBeg[0] = fPadBeg[1] = fCathBeg;
92
93 if (!fgMinuit) fgMinuit = new TMinuit(8);
94
9bbd7f60 95 if (fPlot) fDebug = 1;
c0a16418 96}
97
98//_____________________________________________________________________________
99AliMUONClusterFinderMLEM::~AliMUONClusterFinderMLEM()
100{
101/// Destructor
102 delete fgMinuit; fgMinuit = 0; delete fPixArray; fPixArray = 0;
103// delete fDraw;
104 delete fPreClusterFinder;
c0a16418 105 delete fSplitter;
106 AliInfo(Form("Total clusters %d AddVirtualPad needed %d",
107 fNClusters,fNAddVirtualPads));
108}
109
110//_____________________________________________________________________________
111Bool_t
24935e58 112AliMUONClusterFinderMLEM::Prepare(Int_t detElemId,
113 TClonesArray* pads[2],
114 const AliMpArea& area,
115 const AliMpVSegmentation* seg[2])
c0a16418 116{
117 /// Prepare for clustering
2060b217 118// AliCodeTimerAuto("")
c0a16418 119
120 for ( Int_t i = 0; i < 2; ++i )
121 {
24935e58 122 fSegmentation[i] = seg[i];
c0a16418 123 }
124
125 // Find out the DetElemId
24935e58 126 fDetElemId = detElemId;
c0a16418 127
128 delete fSplitter;
129 fSplitter = new AliMUONClusterSplitterMLEM(fDetElemId,fPixArray);
9bbd7f60 130 fSplitter->SetDebug(fDebug);
c0a16418 131
132 // find out current event number, and reset the cluster number
133 fEventNumber = AliRunLoader::GetRunLoader()->GetEventNumber();
134 fClusterNumber = -1;
135 fClusterList.Delete();
136
b1a19e07 137 AliDebug(3,Form("EVT %d DE %d",fEventNumber,fDetElemId));
c0a16418 138
24935e58 139 if ( fPreClusterFinder->NeedSegmentation() )
140 {
141 return fPreClusterFinder->Prepare(detElemId,pads,area,seg);
142 }
143 else
144 {
145 return fPreClusterFinder->Prepare(detElemId,pads,area);
146 }
c0a16418 147}
148
149//_____________________________________________________________________________
150AliMUONCluster*
151AliMUONClusterFinderMLEM::NextCluster()
152{
153 /// Return next cluster
f6c291ef 154// AliCodeTimerAuto("")
c0a16418 155
c0a16418 156 // if the list of clusters is not void, pick one from there
9bbd7f60 157 TObject* o = fClusterList.At(++fClusterNumber);
158 if ( o != 0x0 ) return static_cast<AliMUONCluster*>(o);
c0a16418 159
160 //FIXME : at this point, must check whether we've used all the digits
161 //from precluster : if not, let the preclustering know about those unused
162 //digits, so it can reuse them
163
164 // if the cluster list is exhausted, we need to go to the next
165 // pre-cluster and treat it
166
167 fPreCluster = fPreClusterFinder->NextCluster();
168
169 if (!fPreCluster)
170 {
171 // we are done
172 return 0x0;
173 }
174
175 fClusterList.Delete(); // reset the list of clusters for this pre-cluster
9bbd7f60 176 fClusterNumber = -1; //AZ
c0a16418 177
178 WorkOnPreCluster();
179
180 // WorkOnPreCluster may have used only part of the pads, so we check that
181 // now, and let the unused pads be reused by the preclustering...
182
2abdae6e 183 Int_t mult = fPreCluster->Multiplicity();
184 for ( Int_t i = 0; i < mult; ++i )
c0a16418 185 {
186 AliMUONPad* pad = fPreCluster->Pad(i);
187 if ( !pad->IsUsed() )
188 {
189 fPreClusterFinder->UsePad(*pad);
190 }
191 }
192
193 return NextCluster();
194}
195
196//_____________________________________________________________________________
197Bool_t
198AliMUONClusterFinderMLEM::WorkOnPreCluster()
199{
200 /// Starting from a precluster, builds a pixel array, and then
201 /// extract clusters from this array
202
2abdae6e 203 // AliCodeTimerAuto("")
9bbd7f60 204
205 if (fDebug) {
206 cout << " *** Event # " << AliRunLoader::GetRunLoader()->GetEventNumber()
207 << " det. elem.: " << fDetElemId << endl;
2abdae6e 208 for (Int_t j = 0; j < fPreCluster->Multiplicity(); ++j) {
9bbd7f60 209 AliMUONPad* pad = fPreCluster->Pad(j);
210 printf(" bbb %3d %1d %8.4f %8.4f %8.4f %8.4f %6.1f %3d %3d %2d %1d %1d \n",
211 j, pad->Cathode(), pad->Coord(0), pad->Coord(1), pad->DX()*2, pad->DY()*2,
212 pad->Charge(), pad->Ix(), pad->Iy(), pad->Status(), pad->IsReal(), pad->IsSaturated());
213 }
214 }
215
c0a16418 216 AliMUONCluster* cluster = CheckPrecluster(*fPreCluster);
c0a16418 217 if (!cluster) return kFALSE;
9bbd7f60 218
c0a16418 219 BuildPixArray(*cluster);
220
221 if ( fPixArray->GetLast() < 0 )
222 {
223 AliDebug(1,"No pixel for the above cluster");
224 delete cluster;
225 return kFALSE;
226 }
227
228 // Use MLEM for cluster finder
229 Int_t nMax = 1, localMax[100], maxPos[100];
230 Double_t maxVal[100];
231
c0a16418 232 Int_t iSimple = 0, nInX = -1, nInY;
233
234 PadsInXandY(*cluster,nInX, nInY);
235
8c718ff8 236 if (nInX < 4 && nInY < 4)
237 {
238 iSimple = 1; // simple cluster
239 }
240 else
c0a16418 241 {
8c718ff8 242 nMax = FindLocalMaxima(fPixArray, localMax, maxVal); // for small clusters just to tag pixels
243 if (nMax > 1) {
244 if (cluster->Multiplicity() <= 50) nMax = 1; // for small clusters
245 if (nMax > 1) TMath::Sort(nMax, maxVal, maxPos, kTRUE); // in descending order
246 }
c0a16418 247 }
248
2abdae6e 249 for (Int_t i = 0; i < nMax; ++i)
c0a16418 250 {
251 if (nMax > 1)
252 {
253 FindCluster(*cluster,localMax, maxPos[i]);
254 }
2abdae6e 255
c0a16418 256 MainLoop(*cluster,iSimple);
2abdae6e 257
c0a16418 258 if (i < nMax-1)
259 {
2abdae6e 260 Int_t mult = cluster->Multiplicity();
261 for (Int_t j = 0; j < mult; ++j)
c0a16418 262 {
263 AliMUONPad* pad = cluster->Pad(j);
8c718ff8 264 //if ( pad->Status() == 0 ) continue; // pad charge was not modified
265 if ( pad->Status() != fgkOver ) continue; // pad was not used
266 //pad->SetStatus(0);
267 pad->SetStatus(fgkZero);
c0a16418 268 pad->RevertCharge(); // use backup charge value
269 }
270 }
271 } // for (Int_t i=0; i<nMax;
8c718ff8 272 //if (nMax > 1) ((TH2D*) gROOT->FindObject("anode"))->Delete();
273 delete ((TH2D*) gROOT->FindObject("anode"));
2abdae6e 274 //TH2D *mlem = (TH2D*) gROOT->FindObject("mlem");
275 //if (mlem) mlem->Delete();
c0a16418 276 delete cluster;
277 return kTRUE;
278}
279
280//_____________________________________________________________________________
281Bool_t
282AliMUONClusterFinderMLEM::Overlap(const AliMUONPad& pad, const AliMUONPad& pixel)
283{
284 /// Check if the pad and the pixel overlaps
285
286 // make a fake pad from the pixel
287 AliMUONPad tmp(pad.DetElemId(),pad.Cathode(),pad.Ix(),pad.Iy(),
288 pixel.Coord(0),pixel.Coord(1),
289 pixel.Size(0),pixel.Size(1),0);
290
291 return AliMUONPad::AreOverlapping(pad,tmp,fgkDecreaseSize);
292}
293
294//_____________________________________________________________________________
295AliMUONCluster*
296AliMUONClusterFinderMLEM::CheckPrecluster(const AliMUONCluster& origCluster)
297{
298 /// Check precluster in order to attempt to simplify it (mostly for
299 /// two-cathode preclusters)
2abdae6e 300
24935e58 301 AliCodeTimerAuto("")
f6c291ef 302
2abdae6e 303 // Disregard small clusters (leftovers from splitting or noise)
304 if ((origCluster.Multiplicity()==1 || origCluster.Multiplicity()==2) &&
305 origCluster.Charge(0)+origCluster.Charge(1) < 10)
c0a16418 306 {
c0a16418 307 return 0x0;
308 }
309
8c718ff8 310 AliMUONCluster* cluster = new AliMUONCluster(origCluster);
c0a16418 311
c0a16418 312 AliDebug(2,"Start of CheckPreCluster=");
2abdae6e 313 //StdoutToAliDebug(2,cluster->Print("full"));
c0a16418 314
c0a16418 315 AliMUONCluster* rv(0x0);
316
2abdae6e 317 if (cluster->Multiplicity(0) && cluster->Multiplicity(1))
c0a16418 318 {
319 rv = CheckPreclusterTwoCathodes(cluster);
320 }
321 else
322 {
f6c291ef 323 rv = cluster;
c0a16418 324 }
c0a16418 325 return rv;
326}
327
c0a16418 328//_____________________________________________________________________________
329AliMUONCluster*
330AliMUONClusterFinderMLEM::CheckPreclusterTwoCathodes(AliMUONCluster* cluster)
331{
71a2d3aa 332 /// Check two-cathode cluster
c0a16418 333
c0a16418 334 Int_t npad = cluster->Multiplicity();
335 Int_t* flags = new Int_t[npad];
2abdae6e 336 for (Int_t j = 0; j < npad; ++j) flags[j] = 0;
c0a16418 337
338 // Check pad overlaps
2abdae6e 339 for ( Int_t i = 0; i < npad; ++i)
c0a16418 340 {
341 AliMUONPad* padi = cluster->Pad(i);
8c718ff8 342 if ( padi->Cathode() != 0 ) continue;
2abdae6e 343 for (Int_t j = i+1; j < npad; ++j)
c0a16418 344 {
345 AliMUONPad* padj = cluster->Pad(j);
8c718ff8 346 if ( padj->Cathode() != 1 ) continue;
c0a16418 347 if ( !AliMUONPad::AreOverlapping(*padi,*padj,fgkDecreaseSize) ) continue;
348 flags[i] = flags[j] = 1; // mark overlapped pads
349 }
350 }
351
352 // Check if all pads overlap
353 Int_t nFlags=0;
2abdae6e 354 for (Int_t i = 0; i < npad; ++i)
c0a16418 355 {
2abdae6e 356 if (!flags[i]) ++nFlags;
c0a16418 357 }
358
359 if (nFlags > 0)
360 {
361 // not all pads overlap.
9bbd7f60 362 if (fDebug) cout << " nFlags: " << nFlags << endl;
363 TObjArray toBeRemoved;
2abdae6e 364 for (Int_t i = 0; i < npad; ++i)
c0a16418 365 {
366 AliMUONPad* pad = cluster->Pad(i);
367 if (flags[i]) continue;
368 Int_t cath = pad->Cathode();
369 Int_t cath1 = TMath::Even(cath);
370 // Check for edge effect (missing pads on the _other_ cathode)
371 AliMpPad mpPad = fSegmentation[cath1]->PadByPosition(pad->Position(),kFALSE);
372 if (!mpPad.IsValid()) continue;
8c718ff8 373 //if (nFlags == 1 && pad->Charge() < fgkZeroSuppression * 3) continue;
374 if (nFlags == 1 && pad->Charge() < 20) continue;
c0a16418 375 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
376 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),pad->Charge()));
9bbd7f60 377 toBeRemoved.AddLast(pad);
c0a16418 378 fPreCluster->Pad(i)->Release();
9bbd7f60 379 }
2abdae6e 380 Int_t nRemove = toBeRemoved.GetEntriesFast();
381 for ( Int_t i = 0; i < nRemove; ++i )
9bbd7f60 382 {
2abdae6e 383 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
c0a16418 384 }
385 }
386
387 // Check correlations of cathode charges
9bbd7f60 388 if ( !cluster->IsSaturated() && cluster->ChargeAsymmetry() > 1 )
c0a16418 389 {
390 // big difference
391 Int_t cathode = cluster->MaxRawChargeCathode();
8c718ff8 392 Int_t imin(-1);
393 Int_t imax(-1);
c0a16418 394 Double_t cmax(0);
395 Double_t cmin(1E9);
396
397 // get min and max pad charges on the cathode opposite to the
398 // max pad (given by MaxRawChargeCathode())
399 //
2abdae6e 400 Int_t mult = cluster->Multiplicity();
401 for ( Int_t i = 0; i < mult; ++i )
c0a16418 402 {
403 AliMUONPad* pad = cluster->Pad(i);
404 if ( pad->Cathode() != cathode || !pad->IsReal() )
405 {
406 // only consider pads in the opposite cathode, and
2abdae6e 407 // only consider real pads (i.e. exclude the virtual ones)
c0a16418 408 continue;
409 }
410 if ( pad->Charge() < cmin )
411 {
412 cmin = pad->Charge();
413 imin = i;
8c718ff8 414 if (imax < 0) {
415 imax = imin;
416 cmax = cmin;
417 }
c0a16418 418 }
8c718ff8 419 else if ( pad->Charge() > cmax )
c0a16418 420 {
421 cmax = pad->Charge();
422 imax = i;
423 }
424 }
425 AliDebug(2,Form("Pad imin,imax %d,%d cmin,cmax %e,%e",
426 imin,imax,cmin,cmax));
427 //
428 // arrange pads according to their distance to the max, normalized
429 // to the pad size
2abdae6e 430 Double_t* dist = new Double_t[mult];
c0a16418 431 Double_t dxMin(1E9);
432 Double_t dyMin(1E9);
433 Double_t dmin(0);
434
435 AliMUONPad* padmax = cluster->Pad(imax);
436
2abdae6e 437 for ( Int_t i = 0; i < mult; ++i )
c0a16418 438 {
439 dist[i] = 0.0;
440 if ( i == imax) continue;
441 AliMUONPad* pad = cluster->Pad(i);
442 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
443 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
444 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
445 dist[i] = TMath::Sqrt(dx*dx+dy*dy);
446 if ( i == imin )
447 {
448 dmin = dist[i] + 1E-3; // distance to the pad with minimum charge
449 dxMin = dx;
450 dyMin = dy;
451 }
452 }
453
2abdae6e 454 TMath::Sort(mult,dist,flags,kFALSE); // in ascending order
455 Double_t xmax(-1), distPrev(999);
c0a16418 456 TObjArray toBeRemoved;
457
2abdae6e 458 for ( Int_t i = 0; i < mult; ++i )
c0a16418 459 {
460 Int_t indx = flags[i];
461 AliMUONPad* pad = cluster->Pad(indx);
462 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
463 if ( dist[indx] > dmin )
464 {
465 // farther than the minimum pad
466 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
467 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
468 dx *= dxMin;
469 dy *= dyMin;
470 if (dx >= 0 && dy >= 0) continue;
471 if (TMath::Abs(dx) > TMath::Abs(dy) && dx >= 0) continue;
472 if (TMath::Abs(dy) > TMath::Abs(dx) && dy >= 0) continue;
473 }
2abdae6e 474 if (dist[indx] > distPrev + 1) break; // overstepping empty pads
c0a16418 475 if ( pad->Charge() <= cmax || TMath::Abs(dist[indx]-xmax) < 1E-3 )
476 {
477 // release pad
478 if (TMath::Abs(dist[indx]-xmax) < 1.e-3)
479 {
480 cmax = TMath::Max(pad->Charge(),cmax);
481 }
482 else
483 {
484 cmax = pad->Charge();
485 }
486 xmax = dist[indx];
2abdae6e 487 distPrev = dist[indx];
c0a16418 488 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
489 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),
490 pad->Charge()));
491
492 toBeRemoved.AddLast(pad);
493 fPreCluster->Pad(indx)->Release();
494 }
495 }
2abdae6e 496 Int_t nRemove = toBeRemoved.GetEntriesFast();
497 for ( Int_t i = 0; i < nRemove; ++i )
c0a16418 498 {
2abdae6e 499 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
c0a16418 500 }
501 delete[] dist;
502 }
503
504 delete[] flags;
505
506 AliDebug(2,"End of CheckPreClusterTwoCathodes=");
2abdae6e 507 //StdoutToAliDebug(2,cluster->Print("full"));
c0a16418 508
509 return cluster;
510}
511
512//_____________________________________________________________________________
513void
514AliMUONClusterFinderMLEM::CheckOverlaps()
515{
516 /// For debug only : check if some pixels overlap...
517
518 Int_t nPix = fPixArray->GetLast()+1;
519 Int_t dummy(0);
520
521 for ( Int_t i = 0; i < nPix; ++i )
522 {
523 AliMUONPad* pixelI = Pixel(i);
524 AliMUONPad pi(dummy,dummy,dummy,dummy,
525 pixelI->Coord(0),pixelI->Coord(1),
526 pixelI->Size(0),pixelI->Size(1),0.0);
527
528 for ( Int_t j = i+1; j < nPix; ++j )
529 {
530 AliMUONPad* pixelJ = Pixel(j);
531 AliMUONPad pj(dummy,dummy,dummy,dummy,
532 pixelJ->Coord(0),pixelJ->Coord(1),
533 pixelJ->Size(0),pixelJ->Size(1),0.0);
534 AliMpArea area;
535
536 if ( AliMUONPad::AreOverlapping(pi,pj,fgkDecreaseSize,area) )
537 {
538 AliInfo(Form("The following 2 pixels (%d and %d) overlap !",i,j));
2abdae6e 539 /*
c0a16418 540 StdoutToAliInfo(pixelI->Print();
541 cout << " Surface = " << pixelI->Size(0)*pixelI->Size(1)*4 << endl;
542 pixelJ->Print();
543 cout << " Surface = " << pixelJ->Size(0)*pixelJ->Size(1)*4 << endl;
544 cout << " Area surface = " << area.Dimensions().X()*area.Dimensions().Y()*4 << endl;
545 cout << "-------" << endl;
546 );
2abdae6e 547 */
c0a16418 548 }
549 }
550 }
551}
552
553//_____________________________________________________________________________
554void AliMUONClusterFinderMLEM::BuildPixArray(AliMUONCluster& cluster)
555{
556 /// Build pixel array for MLEM method
557
558 Int_t npad = cluster.Multiplicity();
559 if (npad<=0)
560 {
561 AliWarning("Got no pad at all ?!");
562 }
563
564 fPixArray->Delete();
8c718ff8 565 BuildPixArrayOneCathode(cluster);
c0a16418 566
567 Int_t nPix = fPixArray->GetLast()+1;
568
2abdae6e 569// AliDebug(2,Form("nPix after BuildPixArray=%d",nPix));
c0a16418 570
2abdae6e 571 if ( nPix > npad )
c0a16418 572 {
2abdae6e 573// AliDebug(2,Form("Will trim number of pixels to number of pads"));
c0a16418 574
575 // Too many pixels - sort and remove pixels with the lowest signal
576 fPixArray->Sort();
8c718ff8 577 for ( Int_t i = npad; i < nPix; ++i )
c0a16418 578 {
579 RemovePixel(i);
580 }
581 fPixArray->Compress();
c0a16418 582 } // if (nPix > npad)
583
584// StdoutToAliDebug(2,cout << "End of BuildPixelArray:" << endl;
585// fPixArray->Print(););
2abdae6e 586 //CheckOverlaps();//FIXME : this is for debug only. Remove it.
c0a16418 587}
588
2abdae6e 589//_____________________________________________________________________________
590void AliMUONClusterFinderMLEM::BuildPixArrayOneCathode(AliMUONCluster& cluster)
591{
592 /// Build the pixel array
593
594// AliDebug(2,Form("cluster.Multiplicity=%d",cluster.Multiplicity()));
595
2abdae6e 596 TVector2 dim = cluster.MinPadDimensions (-1, kFALSE);
597 Double_t width[2] = {dim.X(), dim.Y()}, xy0[2];
8c718ff8 598 Int_t found[2] = {0}, mult = cluster.Multiplicity();
2abdae6e 599
8c718ff8 600 for ( Int_t i = 0; i < mult; ++i) {
2abdae6e 601 AliMUONPad* pad = cluster.Pad(i);
2abdae6e 602 for (Int_t j = 0; j < 2; ++j) {
603 if (found[j] == 0 && TMath::Abs(pad->Size(j)-width[j]) < fgkDistancePrecision) {
604 xy0[j] = pad->Coord(j);
605 found[j] = 1;
606 }
607 }
8c718ff8 608 if (found[0] && found[1]) break;
2abdae6e 609 }
2abdae6e 610
2abdae6e 611 Double_t min[2], max[2];
612 Int_t cath0 = 0, cath1 = 1;
613 if (cluster.Multiplicity(0) == 0) cath0 = 1;
614 else if (cluster.Multiplicity(1) == 0) cath1 = 0;
8c718ff8 615
616 TVector2 leftDown = cluster.Area(cath0).LeftDownCorner();
617 TVector2 rightUp = cluster.Area(cath0).RightUpCorner();
618 min[0] = leftDown.X();
619 min[1] = leftDown.Y();
620 max[0] = rightUp.X();
621 max[1] = rightUp.Y();
622 if (cath1 != cath0) {
623 leftDown = cluster.Area(cath1).LeftDownCorner();
624 rightUp = cluster.Area(cath1).RightUpCorner();
625 min[0] = TMath::Max (min[0], leftDown.X());
626 min[1] = TMath::Max (min[1], leftDown.Y());
627 max[0] = TMath::Min (max[0], rightUp.X());
628 max[1] = TMath::Min (max[1], rightUp.Y());
629 }
2abdae6e 630
631 // Adjust limits
632 //width[0] /= 2; width[1] /= 2; // just for check
633 Int_t nbins[2];
634 for (Int_t i = 0; i < 2; ++i) {
635 Double_t dist = (min[i] - xy0[i]) / width[i] / 2;
8c718ff8 636 if (TMath::Abs(dist) < 1.e-6) dist = -1.e-6;
2abdae6e 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;
643 }
644
645 // Book histogram
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();
650
651 // Fill histogram
2abdae6e 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);
657 }
658
659 // Store pixels
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;
8c718ff8 664 //if (hist2->GetCellContent(i,j) < 1.1 && cluster.Multiplicity(0) &&
665 // cluster.Multiplicity(1)) continue;
666 if (cath0 != cath1) {
667 // Two-sided cluster
668 Double_t cont = hist2->GetCellContent(i,j);
669 if (cont < 999.) continue;
670 if (cont-Int_t(cont/1000.)*1000. < 0.5) continue;
671 }
2abdae6e 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);
676 }
677 }
8c718ff8 678 //*
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);
686 }
687 //*/
2abdae6e 688 //fPixArray->Print();
689 delete hist1;
690 delete hist2;
691}
692
693//_____________________________________________________________________________
694void AliMUONClusterFinderMLEM::PadOverHist(Int_t idir, Int_t ix0, Int_t iy0, AliMUONPad *pad)
695{
696 /// "Span" pad over histogram in the direction idir
697
698 TH2D *hist1 = static_cast<TH2D*> (gROOT->FindObject("Grid"));
699 TH2D *hist2 = static_cast<TH2D*> (gROOT->FindObject("Entries"));
700 TAxis *axis = idir == 0 ? hist1->GetXaxis() : hist1->GetYaxis();
8c718ff8 701 Int_t nbins = axis->GetNbins(), cath = pad->Cathode();
702 Double_t bin = axis->GetBinWidth(1), amask = TMath::Power(1000.,cath*1.);
2abdae6e 703
704 Int_t nbinPad = (Int_t)(pad->Size(idir)/bin*2+fgkDistancePrecision) + 1; // number of bins covered by pad
705
706 for (Int_t i = 0; i < nbinPad; ++i) {
707 Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
708 if (ixy > nbins) break;
709 Double_t lowEdge = axis->GetBinLowEdge(ixy);
710 if (lowEdge + fgkDistancePrecision > pad->Coord(idir) + pad->Size(idir)) break;
711 if (idir == 0) PadOverHist(1, ixy, iy0, pad); // span in the other direction
712 else {
713 // Fill histogram
714 Double_t cont = pad->Charge();
715 if (hist2->GetCellContent(ix0, ixy) > 0.1)
716 cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont);
717 hist1->SetCellContent(ix0, ixy, cont);
8c718ff8 718 //hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+1);
719 hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
2abdae6e 720 }
721 }
722
723 for (Int_t i = -1; i > -nbinPad; --i) {
724 Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
725 if (ixy < 1) break;
726 Double_t upEdge = axis->GetBinUpEdge(ixy);
727 if (upEdge - fgkDistancePrecision < pad->Coord(idir) - pad->Size(idir)) break;
728 if (idir == 0) PadOverHist(1, ixy, iy0, pad); // span in the other direction
729 else {
730 // Fill histogram
731 Double_t cont = pad->Charge();
732 if (hist2->GetCellContent(ix0, ixy) > 0.1)
733 cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont);
734 hist1->SetCellContent(ix0, ixy, cont);
8c718ff8 735 //hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+1);
736 hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
2abdae6e 737 }
738 }
739}
740
c0a16418 741//_____________________________________________________________________________
742void
743AliMUONClusterFinderMLEM::Plot(const char* basename)
744{
745 /// Make a plot and save it as png
746
9bbd7f60 747 return; //AZ
c0a16418 748 if (!fPlot) return;
749
750 TCanvas* c = new TCanvas("MLEM","MLEM",800,600);
751 c->Draw();
752 Draw();
753 c->Modified();
754 c->Update();
755 c->Print(Form("%s.EVT%d.DE%d.CLU%d.png",basename,fEventNumber,
756 fDetElemId,fClusterNumber));
757}
758
759//_____________________________________________________________________________
760void
761AliMUONClusterFinderMLEM::ComputeCoefficients(AliMUONCluster& cluster,
762 Double_t* coef,
763 Double_t* probi)
764{
765 /// Compute coefficients needed for MLEM algorithm
766
767 Int_t nPix = fPixArray->GetLast()+1;
768
2abdae6e 769 //memset(probi,0,nPix*sizeof(Double_t));
770 for (Int_t j = 0; j < nPix; ++j) probi[j] = 0.;
c0a16418 771
2abdae6e 772 Int_t mult = cluster.Multiplicity();
773 for ( Int_t j = 0; j < mult; ++j )
c0a16418 774 {
775 AliMUONPad* pad = cluster.Pad(j);
776 Int_t indx = j*nPix;
777
2abdae6e 778 for ( Int_t ipix = 0; ipix < nPix; ++ipix )
c0a16418 779 {
780 Int_t indx1 = indx + ipix;
8c718ff8 781 //if (pad->Status() < 0)
782 if (pad->Status() != fgkZero)
c0a16418 783 {
784 coef[indx1] = 0;
785 continue;
786 }
787 AliMUONPad* pixPtr = Pixel(ipix);
788 // coef is the charge (given by Mathieson integral) on pad, assuming
789 // the Mathieson is center at pixel.
790 coef[indx1] = fSplitter->ChargeIntegration(pixPtr->Coord(0), pixPtr->Coord(1), *pad);
2abdae6e 791// AliDebug(2,Form("pad=(%d,%d,%e,%e,%e,%e) pix=(%e,%e,%e,%e) coef %e",
792// pad->Ix(),pad->Iy(),
793// pad->X(),pad->Y(),
794// pad->DX(),pad->DY(),
795// pixPtr->Coord(0),pixPtr->Coord(1),
796// pixPtr->Size(0),pixPtr->Size(1),
797// coef[indx1]));
c0a16418 798
799 probi[ipix] += coef[indx1];
800 }
801 }
802}
803
804//_____________________________________________________________________________
805Bool_t AliMUONClusterFinderMLEM::MainLoop(AliMUONCluster& cluster, Int_t iSimple)
806{
807 /// Repeat MLEM algorithm until pixel size becomes sufficiently small
808
2abdae6e 809 // AliCodeTimerAuto("")
f6c291ef 810
c0a16418 811 Int_t nPix = fPixArray->GetLast()+1;
812
813 AliDebug(2,Form("nPix=%d iSimple=%d, precluster=",nPix,iSimple));
2abdae6e 814 //StdoutToAliDebug(2,cluster.Print("full"););
c0a16418 815
816 if ( nPix < 0 )
817 {
818 AliDebug(1,"No pixels, why am I here then ?");
819 }
820
821 AddVirtualPad(cluster); // add virtual pads if necessary
822
823 Int_t npadTot = cluster.Multiplicity();
824 Int_t npadOK = 0;
825 for (Int_t i = 0; i < npadTot; ++i)
826 {
8c718ff8 827 //if (cluster.Pad(i)->Status() == 0) ++npadOK;
828 if (cluster.Pad(i)->Status() == fgkZero) ++npadOK;
c0a16418 829 }
830
831 TH2D* mlem(0x0);
832 Double_t* coef(0x0);
833 Double_t* probi(0x0);
8c718ff8 834 Int_t lc(0); // loop counter
c0a16418 835
2abdae6e 836 //Plot("mlem.start");
8c718ff8 837 AliMUONPad* pixPtr = Pixel(0);
838 Double_t xylim[4] = {pixPtr->X(), -pixPtr->X(), pixPtr->Y(), -pixPtr->Y()};
839
c0a16418 840 while (1)
841 {
842 ++lc;
843 mlem = (TH2D*) gROOT->FindObject("mlem");
844 delete mlem;
845
846 AliDebug(2,Form("lc %d nPix %d(%d) npadTot %d npadOK %d",lc,nPix,fPixArray->GetLast()+1,npadTot,npadOK));
847 AliDebug(2,Form("EVT%d PixArray=",fEventNumber));
2abdae6e 848 //StdoutToAliDebug(2,fPixArray->Print("","full"));
c0a16418 849
850 coef = new Double_t [npadTot*nPix];
851 probi = new Double_t [nPix];
852
853 // Calculate coefficients and pixel visibilities
854 ComputeCoefficients(cluster,coef,probi);
855
2abdae6e 856 for (Int_t ipix = 0; ipix < nPix; ++ipix)
c0a16418 857 {
858 if (probi[ipix] < 0.01)
859 {
860 AliMUONPad* pixel = Pixel(ipix);
861 AliDebug(2,Form("Setting the following pixel to invisible as its probi<0.01:"));
2abdae6e 862 //StdoutToAliDebug(2,cout << Form(" -- ipix %3d --- "); pixel->Print(););
c0a16418 863 pixel->SetCharge(0); // "invisible" pixel
864 }
865 }
866
867 // MLEM algorithm
868 Mlem(cluster,coef, probi, 15);
869
8c718ff8 870 // Find histogram limits for the 1'st pass only - for others computed below
871 if (lc == 1) {
872 for ( Int_t ipix = 1; ipix < nPix; ++ipix )
873 {
874 pixPtr = Pixel(ipix);
875 for ( Int_t i = 0; i < 2; ++i )
876 {
877 Int_t indx = i * 2;
878 if (pixPtr->Coord(i) < xylim[indx]) xylim[indx] = pixPtr->Coord(i);
879 else if (-pixPtr->Coord(i) < xylim[indx+1]) xylim[indx+1] = -pixPtr->Coord(i);
880 }
881 }
882 } else pixPtr = Pixel(0);
883
2abdae6e 884 for (Int_t i = 0; i < 4; i++)
c0a16418 885 {
886 xylim[i] -= pixPtr->Size(i/2);
887 }
c0a16418 888
889 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
890 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
891
2abdae6e 892 //StdoutToAliDebug(2,cout << "pixel used for nx,ny computation : "; pixPtr->Print(););
c0a16418 893 AliDebug(2,Form("lc %d pixPtr size = %e,%e nx,ny=%d,%d xylim=%e,%e,%e,%e",
894 lc,pixPtr->Size(0),pixPtr->Size(1),nx,ny,
895 xylim[0],-xylim[1],xylim[2],-xylim[3]
896 ));
897
898 mlem = new TH2D("mlem","mlem",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
899
2abdae6e 900 for (Int_t ipix = 0; ipix < nPix; ++ipix)
c0a16418 901 {
902 AliMUONPad* pixPtr = Pixel(ipix);
903 mlem->Fill(pixPtr->Coord(0),pixPtr->Coord(1),pixPtr->Charge());
904 }
905
906 // Check if the total charge of pixels is too low
907 Double_t qTot = 0;
2abdae6e 908 for ( Int_t i = 0; i < nPix; ++i)
c0a16418 909 {
910 qTot += Pixel(i)->Charge();
911 }
912
913 if ( qTot < 1.e-4 || ( npadOK < 3 && qTot < 7 ) )
914 {
915 AliDebug(1,Form("Deleting the above cluster (charge %e too low, npadOK=%d)",qTot,npadOK));
916 delete [] coef;
917 delete [] probi;
c0a16418 918 fPixArray->Delete();
2abdae6e 919 for ( Int_t i = 0; i < npadTot; ++i)
c0a16418 920 {
921 AliMUONPad* pad = cluster.Pad(i);
8c718ff8 922 //if ( pad->Status() == 0) pad->SetStatus(-1);
923 if ( pad->Status() == fgkZero) pad->SetStatus(fgkOver);
c0a16418 924 }
925 return kFALSE;
926 }
927
928 if (iSimple)
929 {
930 // Simple cluster - skip further passes thru EM-procedure
931 Simple(cluster);
932 delete [] coef;
933 delete [] probi;
c0a16418 934 fPixArray->Delete();
935 return kTRUE;
936 }
937
938 // Calculate position of the center-of-gravity around the maximum pixel
939 Double_t xyCOG[2];
940 FindCOG(mlem, xyCOG);
941
942 if (TMath::Min(pixPtr->Size(0),pixPtr->Size(1)) < 0.07 &&
943 pixPtr->Size(0) > pixPtr->Size(1)) break;
944
945 // Sort pixels according to the charge
8c718ff8 946 MaskPeaks(1); // mask local maxima
c0a16418 947 fPixArray->Sort();
8c718ff8 948 MaskPeaks(0); // unmask local maxima
c0a16418 949 Double_t pixMin = 0.01*Pixel(0)->Charge();
950 pixMin = TMath::Min(pixMin,50.);
951
952 // Decrease pixel size and shift pixels to make them centered at
953 // the maximum one
954 Int_t indx = (pixPtr->Size(0)>pixPtr->Size(1)) ? 0 : 1;
955 Int_t ix(1);
956 Double_t width = 0;
957 Double_t shift[2] = { 0.0, 0.0 };
2abdae6e 958 for (Int_t i = 0; i < 4; ++i) xylim[i] = 999;
c0a16418 959 Int_t nPix1 = nPix;
960 nPix = 0;
2abdae6e 961 for (Int_t ipix = 0; ipix < nPix1; ++ipix)
c0a16418 962 {
963 AliMUONPad* pixPtr = Pixel(ipix);
964 if ( nPix >= npadOK // too many pixels already
965 ||
8c718ff8 966 pixPtr->Charge() < pixMin && pixPtr->Status() != fgkMustKeep // too low charge
c0a16418 967 )
968 {
969 RemovePixel(ipix);
970 continue;
971 }
2abdae6e 972 for (Int_t i = 0; i < 2; ++i)
c0a16418 973 {
974 if (!i)
975 {
976 pixPtr->SetCharge(10);
977 pixPtr->SetSize(indx, pixPtr->Size(indx)/2);
978 width = -pixPtr->Size(indx);
979 pixPtr->Shift(indx, width);
980 // Shift pixel position
981 if (ix)
982 {
983 ix = 0;
2abdae6e 984 for (Int_t j = 0; j < 2; ++j)
c0a16418 985 {
986 shift[j] = pixPtr->Coord(j) - xyCOG[j];
987 shift[j] -= ((Int_t)(shift[j]/pixPtr->Size(j)/2))*pixPtr->Size(j)*2;
988 }
989 } // if (ix)
990 pixPtr->Shift(0, -shift[0]);
991 pixPtr->Shift(1, -shift[1]);
8c718ff8 992 ++nPix;
c0a16418 993 }
8c718ff8 994 else if (nPix < npadOK)
c0a16418 995 {
996 pixPtr = new AliMUONPad(*pixPtr);
997 pixPtr->Shift(indx, -2*width);
8c718ff8 998 pixPtr->SetStatus(fgkZero);
c0a16418 999 fPixArray->Add(pixPtr);
8c718ff8 1000 ++nPix;
c0a16418 1001 }
8c718ff8 1002 else continue; // skip adjustment of histo limits
1003 for (Int_t j = 0; j < 4; ++j)
c0a16418 1004 {
8c718ff8 1005 xylim[j] = TMath::Min (xylim[j], (j%2 ? -1 : 1)*pixPtr->Coord(j/2));
c0a16418 1006 }
1007 } // for (Int_t i=0; i<2;
c0a16418 1008 } // for (Int_t ipix=0;
1009
1010 fPixArray->Compress();
c0a16418 1011
1012 AliDebug(2,Form("After shift:"));
2abdae6e 1013 //StdoutToAliDebug(2,fPixArray->Print("","full"););
1014 //Plot(Form("mlem.lc%d",lc+1));
c0a16418 1015
1016 AliDebug(2,Form(" xyCOG=%9.6f %9.6f xylim=%9.6f,%9.6f,%9.6f,%9.6f",
1017 xyCOG[0],xyCOG[1],
1018 xylim[0],xylim[1],
1019 xylim[2],xylim[3]));
1020
8c718ff8 1021 if (nPix < npadOK)
c0a16418 1022 {
1023 AliMUONPad* pixPtr = Pixel(0);
8c718ff8 1024 // add pixels if the maximum is at the limit of pixel area:
c0a16418 1025 // start from Y-direction
1026 Int_t j = 0;
2abdae6e 1027 for (Int_t i = 3; i > -1; --i)
c0a16418 1028 {
1029 if (nPix < npadOK &&
1030 TMath::Abs((i%2 ? -1 : 1)*xylim[i]-xyCOG[i/2]) < pixPtr->Size(i/2))
1031 {
8c718ff8 1032 //AliMUONPad* p = static_cast<AliMUONPad*>(pixPtr->Clone());
1033 AliMUONPad* p = new AliMUONPad(*pixPtr);
c0a16418 1034 p->SetCoord(i/2, xyCOG[i/2]+(i%2 ? 2:-2)*pixPtr->Size(i/2));
8c718ff8 1035 xylim[i] = p->Coord(i/2) * (i%2 ? -1 : 1); // update histo limits
c0a16418 1036 j = TMath::Even (i/2);
1037 p->SetCoord(j, xyCOG[j]);
1038 AliDebug(2,Form("Adding pixel on the edge (i=%d) ",i));
2abdae6e 1039 //StdoutToAliDebug(2,cout << " ---- ";
1040 // p->Print("corners"););
c0a16418 1041 fPixArray->Add(p);
1042 ++nPix;
1043 }
1044 }
1045 }
c0a16418 1046 nPix = fPixArray->GetEntriesFast();
1047 delete [] coef;
1048 delete [] probi;
c0a16418 1049 } // while (1)
1050
1051 AliDebug(2,Form("At the end of while loop nPix=%d : ",fPixArray->GetLast()+1));
2abdae6e 1052 //StdoutToAliDebug(2,fPixArray->Print("","full"););
c0a16418 1053
1054 // remove pixels with low signal or low visibility
1055 // Cuts are empirical !!!
1056 Double_t thresh = TMath::Max (mlem->GetMaximum()/100.,1.);
1057 thresh = TMath::Min (thresh,50.);
c0a16418 1058 Double_t charge = 0;
1059
c0a16418 1060 // Mark pixels which should be removed
2abdae6e 1061 for (Int_t i = 0; i < nPix; ++i)
c0a16418 1062 {
1063 AliMUONPad* pixPtr = Pixel(i);
1064 charge = pixPtr->Charge();
1065 if (charge < thresh)
1066 {
1067 pixPtr->SetCharge(-charge);
1068 }
1069 }
1070
1071 // Move charge of removed pixels to their nearest neighbour (to keep total charge the same)
1072 Int_t near = 0;
2abdae6e 1073 for (Int_t i = 0; i < nPix; ++i)
c0a16418 1074 {
1075 AliMUONPad* pixPtr = Pixel(i);
1076 charge = pixPtr->Charge();
1077 if (charge > 0) continue;
1078 near = FindNearest(pixPtr);
1079 pixPtr->SetCharge(0);
1080 probi[i] = 0; // make it "invisible"
1081 AliMUONPad* pnear = Pixel(near);
1082 pnear->SetCharge(pnear->Charge() + (-charge));
1083 }
1084 Mlem(cluster,coef,probi,2);
1085
1086 AliDebug(2,Form("Before splitting nPix=%d EVT %d DE %d",fPixArray->GetLast()+1,fEventNumber,fDetElemId));
2abdae6e 1087 //StdoutToAliDebug(2,fPixArray->Print("","full"););
1088 //Plot("mlem.beforesplit");
c0a16418 1089
8c718ff8 1090 // Update histogram
2abdae6e 1091 for (Int_t i = 0; i < nPix; ++i)
c0a16418 1092 {
1093 AliMUONPad* pixPtr = Pixel(i);
1094 Int_t ix = mlem->GetXaxis()->FindBin(pixPtr->Coord(0));
1095 Int_t iy = mlem->GetYaxis()->FindBin(pixPtr->Coord(1));
1096 mlem->SetBinContent(ix, iy, pixPtr->Charge());
1097 }
1098
1099 // Try to split into clusters
1100 Bool_t ok = kTRUE;
1101 if (mlem->GetSum() < 1)
1102 {
1103 ok = kFALSE;
1104 }
1105 else
1106 {
1107 fSplitter->Split(cluster,mlem,coef,fClusterList);
1108 }
1109
1110 delete [] coef;
1111 delete [] probi;
c0a16418 1112 fPixArray->Delete();
1113
1114 return ok;
1115}
1116
8c718ff8 1117//_____________________________________________________________________________
1118void AliMUONClusterFinderMLEM::MaskPeaks(Int_t mask)
1119{
1120 /// Mask/unmask pixels corresponding to local maxima (add/subtract 10000 to their charge
1121 /// - to avoid loosing low charge pixels after sorting)
1122
1123 for (Int_t i = 0; i < fPixArray->GetEntriesFast(); ++i) {
1124 AliMUONPad* pix = Pixel(i);
1125 if (pix->Status() == fgkMustKeep) {
1126 if (mask == 1) pix->SetCharge(pix->Charge()+10000.);
1127 else pix->SetCharge(pix->Charge()-10000.);
1128 }
1129 }
1130}
1131
c0a16418 1132//_____________________________________________________________________________
1133void AliMUONClusterFinderMLEM::Mlem(AliMUONCluster& cluster,
1134 Double_t* coef, Double_t* probi,
1135 Int_t nIter)
1136{
1137 /// Use MLEM to find pixel charges
1138
1139 Int_t nPix = fPixArray->GetEntriesFast();
1140
1141 Int_t npad = cluster.Multiplicity();
1142
1143 Double_t* probi1 = new Double_t[nPix];
2abdae6e 1144 Double_t probMax = TMath::MaxElement(nPix,probi);
c0a16418 1145
2abdae6e 1146 for (Int_t iter = 0; iter < nIter; ++iter)
c0a16418 1147 {
1148 // Do iterations
2abdae6e 1149 for (Int_t ipix = 0; ipix < nPix; ++ipix)
c0a16418 1150 {
9bbd7f60 1151 Pixel(ipix)->SetChargeBackup(0);
c0a16418 1152 // Correct each pixel
1153 if (probi[ipix] < 0.01) continue; // skip "invisible" pixel
1154 Double_t sum = 0;
1155 probi1[ipix] = probMax;
2abdae6e 1156 for (Int_t j = 0; j < npad; ++j)
c0a16418 1157 {
1158 AliMUONPad* pad = cluster.Pad(j);
8c718ff8 1159 //if (pad->Status() < 0) continue;
1160 if (pad->Status() != fgkZero) continue;
c0a16418 1161 Double_t sum1 = 0;
1162 Int_t indx1 = j*nPix;
1163 Int_t indx = indx1 + ipix;
1164 // Calculate expectation
2abdae6e 1165 for (Int_t i = 0; i < nPix; ++i)
c0a16418 1166 {
1167 sum1 += Pixel(i)->Charge()*coef[indx1+i];
2abdae6e 1168 //cout << i << " " << Pixel(i)->Charge() << " " << coef[indx1+i] << endl;
c0a16418 1169 }
2abdae6e 1170 if ( pad->IsSaturated() && sum1 > pad->Charge() )
c0a16418 1171 {
c0a16418 1172 // correct for pad charge overflows
1173 probi1[ipix] -= coef[indx];
1174 continue;
2abdae6e 1175 //sum1 = pad->Charge();
c0a16418 1176 }
1177
2abdae6e 1178 if (sum1 > 1.e-6) sum += pad->Charge()*coef[indx]/sum1;
1179 //if (coef[indx] > 1.e-6) sum += pad->Charge()*coef[indx]/sum1;
c0a16418 1180 } // for (Int_t j=0;
1181 AliMUONPad* pixPtr = Pixel(ipix);
1182 if (probi1[ipix] > 1.e-6)
1183 {
9bbd7f60 1184 //AZ pixPtr->SetCharge(pixPtr->Charge()*sum/probi1[ipix]);
1185 pixPtr->SetChargeBackup(pixPtr->Charge()*sum/probi1[ipix]);
c0a16418 1186 }
2abdae6e 1187 //cout << " xxx " << ipix << " " << pixPtr->Charge() << " " << pixPtr->ChargeBackup() << " " << sum << " " << probi1[ipix] << endl;
c0a16418 1188 } // for (Int_t ipix=0;
2abdae6e 1189 Double_t qTot = 0;
9bbd7f60 1190 for (Int_t i = 0; i < nPix; ++i) {
1191 AliMUONPad* pixPtr = Pixel(i);
2abdae6e 1192 pixPtr->RevertCharge();
1193 qTot += pixPtr->Charge();
1194 }
1195 if (qTot < 1.e-6) {
1196 // Can happen in clusters with large number of overflows - speeding up
1197 delete [] probi1;
1198 return;
9bbd7f60 1199 }
c0a16418 1200 } // for (Int_t iter=0;
1201 delete [] probi1;
1202}
1203
1204//_____________________________________________________________________________
1205void AliMUONClusterFinderMLEM::FindCOG(TH2D *mlem, Double_t *xyc)
1206{
1207 /// Calculate position of the center-of-gravity around the maximum pixel
1208
1209 Int_t ixmax, iymax, ix, nsumx=0, nsumy=0, nsum=0;
1210 Int_t i1 = -9, j1 = -9;
1211 mlem->GetMaximumBin(ixmax,iymax,ix);
1212 Int_t nx = mlem->GetNbinsX();
1213 Int_t ny = mlem->GetNbinsY();
1214 Double_t thresh = mlem->GetMaximum()/10;
1215 Double_t x, y, cont, xq=0, yq=0, qq=0;
1216
2abdae6e 1217 Int_t ie = TMath::Min(ny,iymax+1), je = TMath::Min(nx,ixmax+1);
1218 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
c0a16418 1219 y = mlem->GetYaxis()->GetBinCenter(i);
2abdae6e 1220 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
c0a16418 1221 cont = mlem->GetCellContent(j,i);
1222 if (cont < thresh) continue;
1223 if (i != i1) {i1 = i; nsumy++;}
1224 if (j != j1) {j1 = j; nsumx++;}
1225 x = mlem->GetXaxis()->GetBinCenter(j);
1226 xq += x*cont;
1227 yq += y*cont;
1228 qq += cont;
1229 nsum++;
1230 }
1231 }
1232
1233 Double_t cmax = 0;
1234 Int_t i2 = 0, j2 = 0;
1235 x = y = 0;
1236 if (nsumy == 1) {
1237 // one bin in Y - add one more (with the largest signal)
2abdae6e 1238 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
c0a16418 1239 if (i == iymax) continue;
2abdae6e 1240 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
c0a16418 1241 cont = mlem->GetCellContent(j,i);
1242 if (cont > cmax) {
1243 cmax = cont;
1244 x = mlem->GetXaxis()->GetBinCenter(j);
1245 y = mlem->GetYaxis()->GetBinCenter(i);
1246 i2 = i;
1247 j2 = j;
1248 }
1249 }
1250 }
1251 xq += x*cmax;
1252 yq += y*cmax;
1253 qq += cmax;
1254 if (i2 != i1) nsumy++;
1255 if (j2 != j1) nsumx++;
1256 nsum++;
1257 } // if (nsumy == 1)
1258
1259 if (nsumx == 1) {
1260 // one bin in X - add one more (with the largest signal)
1261 cmax = x = y = 0;
2abdae6e 1262 for (Int_t j = TMath::Max(1,ixmax-1); j <= je; ++j) {
c0a16418 1263 if (j == ixmax) continue;
2abdae6e 1264 for (Int_t i = TMath::Max(1,iymax-1); i <= ie; ++i) {
c0a16418 1265 cont = mlem->GetCellContent(j,i);
1266 if (cont > cmax) {
1267 cmax = cont;
1268 x = mlem->GetXaxis()->GetBinCenter(j);
1269 y = mlem->GetYaxis()->GetBinCenter(i);
1270 i2 = i;
1271 j2 = j;
1272 }
1273 }
1274 }
1275 xq += x*cmax;
1276 yq += y*cmax;
1277 qq += cmax;
1278 if (i2 != i1) nsumy++;
1279 if (j2 != j1) nsumx++;
1280 nsum++;
1281 } // if (nsumx == 1)
1282
1283 xyc[0] = xq/qq; xyc[1] = yq/qq;
1284 AliDebug(2,Form("x,y COG = %e,%e",xyc[0],xyc[1]));
1285}
1286
1287//_____________________________________________________________________________
1288Int_t AliMUONClusterFinderMLEM::FindNearest(AliMUONPad *pixPtr0)
1289{
1290/// Find the pixel nearest to the given one
1291/// (algorithm may be not very efficient)
1292
1293 Int_t nPix = fPixArray->GetEntriesFast(), imin = 0;
1294 Double_t rmin = 99999, dx = 0, dy = 0, r = 0;
1295 Double_t xc = pixPtr0->Coord(0), yc = pixPtr0->Coord(1);
1296 AliMUONPad *pixPtr;
1297
2abdae6e 1298 for (Int_t i = 0; i < nPix; ++i) {
c0a16418 1299 pixPtr = (AliMUONPad*) fPixArray->UncheckedAt(i);
8c718ff8 1300 if (pixPtr == pixPtr0 || pixPtr->Charge() < 0.5) continue;
c0a16418 1301 dx = (xc - pixPtr->Coord(0)) / pixPtr->Size(0);
1302 dy = (yc - pixPtr->Coord(1)) / pixPtr->Size(1);
1303 r = dx *dx + dy * dy;
1304 if (r < rmin) { rmin = r; imin = i; }
1305 }
1306 return imin;
1307}
1308
c0a16418 1309//_____________________________________________________________________________
1310void
1311AliMUONClusterFinderMLEM::Paint(Option_t*)
1312{
1313 /// Paint cluster and pixels
1314
1315 AliMpArea area(fPreCluster->Area());
1316
1317 gPad->Range(area.LeftBorder(),area.DownBorder(),area.RightBorder(),area.UpBorder());
1318
1319 gVirtualX->SetFillStyle(1001);
1320 gVirtualX->SetFillColor(3);
1321 gVirtualX->SetLineColor(3);
1322
1323 Double_t s(1.0);
1324
1325 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1326 {
1327 AliMUONPad* pixel = Pixel(i);
1328
1329 gPad->PaintBox(pixel->Coord(0)-pixel->Size(0)*s,
1330 pixel->Coord(1)-pixel->Size(1)*s,
1331 pixel->Coord(0)+pixel->Size(0)*s,
1332 pixel->Coord(1)+pixel->Size(1)*s);
1333
1334// for ( Int_t sign = -1; sign < 2; sign +=2 )
1335// {
1336// gPad->PaintLine(pixel->Coord(0) - pixel->Size(0),
1337// pixel->Coord(1) + sign*pixel->Size(1),
1338// pixel->Coord(0) + pixel->Size(0),
1339// pixel->Coord(1) - sign*pixel->Size(1)
1340// );
1341// }
1342 }
1343
1344
1345 gVirtualX->SetFillStyle(0);
1346
1347 fPreCluster->Paint();
1348
1349 gVirtualX->SetLineColor(1);
1350 gVirtualX->SetLineWidth(2);
1351 gVirtualX->SetFillStyle(0);
1352 gVirtualX->SetTextColor(1);
1353 gVirtualX->SetTextAlign(22);
1354
1355 for ( Int_t i = 0; i <= fPixArray->GetLast(); ++i)
1356 {
1357 AliMUONPad* pixel = Pixel(i);
1358 gPad->PaintBox(pixel->Coord(0)-pixel->Size(0),
1359 pixel->Coord(1)-pixel->Size(1),
1360 pixel->Coord(0)+pixel->Size(0),
1361 pixel->Coord(1)+pixel->Size(1));
1362 gVirtualX->SetTextSize(pixel->Size(0)*60);
1363
1364 gPad->PaintText(pixel->Coord(0),pixel->Coord(1),Form("%d",(Int_t)(pixel->Charge())));
1365 }
1366}
1367
1368//_____________________________________________________________________________
1369Int_t AliMUONClusterFinderMLEM::FindLocalMaxima(TObjArray *pixArray, Int_t *localMax, Double_t *maxVal)
1370{
1371/// Find local maxima in pixel space for large preclusters in order to
1372/// try to split them into smaller pieces (to speed up the MLEM procedure)
1373/// or to find additional fitting seeds if clusters were not completely resolved
1374
1375 AliDebug(1,Form("nPix=%d",pixArray->GetLast()+1));
1376
1377 TH2D *hist = NULL;
1378 //if (pixArray == fPixArray) hist = (TH2D*) gROOT->FindObject("anode");
1379 //else { hist = (TH2D*) gROOT->FindObject("anode1"); cout << hist << endl; }
1380 //if (hist) hist->Delete();
1381
1382 Double_t xylim[4] = {999, 999, 999, 999};
1383
1384 Int_t nPix = pixArray->GetEntriesFast();
1385 AliMUONPad *pixPtr = 0;
2abdae6e 1386 for (Int_t ipix = 0; ipix < nPix; ++ipix) {
c0a16418 1387 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
2abdae6e 1388 for (Int_t i = 0; i < 4; ++i)
c0a16418 1389 xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
1390 }
2abdae6e 1391 for (Int_t i = 0; i < 4; ++i) xylim[i] -= pixPtr->Size(i/2);
c0a16418 1392
1393 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
1394 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
1395 if (pixArray == fPixArray) hist = new TH2D("anode","anode",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1396 else hist = new TH2D("anode1","anode1",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
2abdae6e 1397 for (Int_t ipix = 0; ipix < nPix; ++ipix) {
c0a16418 1398 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
1399 hist->Fill(pixPtr->Coord(0), pixPtr->Coord(1), pixPtr->Charge());
1400 }
1401// if (fDraw && pixArray == fPixArray) fDraw->DrawHist("c2", hist);
1402
2abdae6e 1403 Int_t nMax = 0, indx, nxy = ny * nx;
1404 Int_t *isLocalMax = new Int_t[nxy];
1405 for (Int_t i = 0; i < nxy; ++i) isLocalMax[i] = 0;
c0a16418 1406
2abdae6e 1407 for (Int_t i = 1; i <= ny; ++i) {
c0a16418 1408 indx = (i-1) * nx;
2abdae6e 1409 for (Int_t j = 1; j <= nx; ++j) {
c0a16418 1410 if (hist->GetCellContent(j,i) < 0.5) continue;
1411 //if (isLocalMax[indx+j-1] < 0) continue;
1412 if (isLocalMax[indx+j-1] != 0) continue;
1413 FlagLocalMax(hist, i, j, isLocalMax);
1414 }
1415 }
1416
2abdae6e 1417 for (Int_t i = 1; i <= ny; ++i) {
c0a16418 1418 indx = (i-1) * nx;
2abdae6e 1419 for (Int_t j = 1; j <= nx; ++j) {
c0a16418 1420 if (isLocalMax[indx+j-1] > 0) {
1421 localMax[nMax] = indx + j - 1;
1422 maxVal[nMax++] = hist->GetCellContent(j,i);
8c718ff8 1423 ((AliMUONPad*)fSplitter->BinToPix(hist, j, i))->SetStatus(fgkMustKeep);
c0a16418 1424 if (nMax > 99) AliFatal(" Too many local maxima !!!");
1425 }
1426 }
1427 }
1428 if (fDebug) cout << " Local max: " << nMax << endl;
2abdae6e 1429 delete [] isLocalMax;
c0a16418 1430 return nMax;
1431}
1432
1433//_____________________________________________________________________________
1434void AliMUONClusterFinderMLEM::FlagLocalMax(TH2D *hist, Int_t i, Int_t j, Int_t *isLocalMax)
1435{
1436/// Flag pixels (whether or not local maxima)
1437
1438 Int_t nx = hist->GetNbinsX();
1439 Int_t ny = hist->GetNbinsY();
1440 Int_t cont = TMath::Nint (hist->GetCellContent(j,i));
1441 Int_t cont1 = 0, indx = (i-1)*nx+j-1, indx1 = 0, indx2 = 0;
1442
2abdae6e 1443 Int_t ie = i + 2, je = j + 2;
1444 for (Int_t i1 = i-1; i1 < ie; ++i1) {
c0a16418 1445 if (i1 < 1 || i1 > ny) continue;
1446 indx1 = (i1 - 1) * nx;
2abdae6e 1447 for (Int_t j1 = j-1; j1 < je; ++j1) {
c0a16418 1448 if (j1 < 1 || j1 > nx) continue;
1449 if (i == i1 && j == j1) continue;
1450 indx2 = indx1 + j1 - 1;
1451 cont1 = TMath::Nint (hist->GetCellContent(j1,i1));
1452 if (cont < cont1) { isLocalMax[indx] = -1; return; }
1453 else if (cont > cont1) isLocalMax[indx2] = -1;
1454 else { // the same charge
1455 isLocalMax[indx] = 1;
1456 if (isLocalMax[indx2] == 0) {
1457 FlagLocalMax(hist, i1, j1, isLocalMax);
1458 if (isLocalMax[indx2] < 0) { isLocalMax[indx] = -1; return; }
1459 else isLocalMax[indx2] = -1;
1460 }
1461 }
1462 }
1463 }
1464 isLocalMax[indx] = 1; // local maximum
1465}
1466
1467//_____________________________________________________________________________
1468void AliMUONClusterFinderMLEM::FindCluster(AliMUONCluster& cluster,
1469 Int_t *localMax, Int_t iMax)
1470{
1471/// Find pixel cluster around local maximum \a iMax and pick up pads
1472/// overlapping with it
1473
1474 TH2D *hist = (TH2D*) gROOT->FindObject("anode");
2abdae6e 1475 /* Just for check
1476 TCanvas* c = new TCanvas("Anode","Anode",800,600);
1477 c->cd();
1478 hist->Draw("lego1Fb"); // debug
1479 c->Update();
1480 Int_t tmp;
1481 cin >> tmp;
1482 */
c0a16418 1483 Int_t nx = hist->GetNbinsX();
1484 Int_t ny = hist->GetNbinsY();
1485 Int_t ic = localMax[iMax] / nx + 1;
1486 Int_t jc = localMax[iMax] % nx + 1;
2abdae6e 1487 Int_t nxy = ny * nx;
1488 Bool_t *used = new Bool_t[nxy];
1489 for (Int_t i = 0; i < nxy; ++i) used[i] = kFALSE;
c0a16418 1490
1491 // Drop all pixels from the array - pick up only the ones from the cluster
1492 fPixArray->Delete();
1493
1494 Double_t wx = hist->GetXaxis()->GetBinWidth(1)/2;
1495 Double_t wy = hist->GetYaxis()->GetBinWidth(1)/2;
1496 Double_t yc = hist->GetYaxis()->GetBinCenter(ic);
1497 Double_t xc = hist->GetXaxis()->GetBinCenter(jc);
1498 Double_t cont = hist->GetCellContent(jc,ic);
1499 fPixArray->Add(new AliMUONPad (xc, yc, wx, wy, cont));
1500 used[(ic-1)*nx+jc-1] = kTRUE;
2abdae6e 1501 AddBinSimple(hist, ic, jc);
1502 //fSplitter->AddBin(hist, ic, jc, 1, used, (TObjArray*)0); // recursive call
c0a16418 1503
1504 Int_t nPix = fPixArray->GetEntriesFast();
1505 Int_t npad = cluster.Multiplicity();
1506
2abdae6e 1507 for (Int_t i = 0; i < nPix; ++i)
c0a16418 1508 {
1509 AliMUONPad* pixPtr = Pixel(i);
1510 pixPtr->SetSize(0,wx);
1511 pixPtr->SetSize(1,wy);
1512 }
1513
1514 // Pick up pads which overlap with found pixels
2abdae6e 1515 for (Int_t i = 0; i < npad; ++i)
c0a16418 1516 {
8c718ff8 1517 //cluster.Pad(i)->SetStatus(-1);
1518 cluster.Pad(i)->SetStatus(fgkOver); // just the dirty trick
c0a16418 1519 }
1520
2abdae6e 1521 for (Int_t i = 0; i < nPix; ++i)
c0a16418 1522 {
1523 AliMUONPad* pixPtr = Pixel(i);
2abdae6e 1524 for (Int_t j = 0; j < npad; ++j)
c0a16418 1525 {
1526 AliMUONPad* pad = cluster.Pad(j);
8c718ff8 1527 //if (pad->Status() == 0) continue;
1528 if (pad->Status() == fgkZero) continue;
c0a16418 1529 if ( Overlap(*pad,*pixPtr) )
1530 {
8c718ff8 1531 //pad->SetStatus(0);
1532 pad->SetStatus(fgkZero);
2abdae6e 1533 if (fDebug) { cout << j << " "; pad->Print("full"); }
c0a16418 1534 }
1535 }
1536 }
1537
2abdae6e 1538 delete [] used;
1539}
1540
1541//_____________________________________________________________________________
1542void
1543AliMUONClusterFinderMLEM::AddBinSimple(TH2D *mlem, Int_t ic, Int_t jc)
1544{
1545 /// Add adjacent bins (+-1 in X and Y) to the cluster
1546
1547 Int_t nx = mlem->GetNbinsX();
1548 Int_t ny = mlem->GetNbinsY();
1549 Double_t cont1, cont = mlem->GetCellContent(jc,ic);
1550 AliMUONPad *pixPtr = 0;
1551
1552 Int_t ie = TMath::Min(ic+1,ny), je = TMath::Min(jc+1,nx);
1553 for (Int_t i = TMath::Max(ic-1,1); i <= ie; ++i) {
1554 for (Int_t j = TMath::Max(jc-1,1); j <= je; ++j) {
1555 cont1 = mlem->GetCellContent(j,i);
1556 if (cont1 > cont) continue;
1557 if (cont1 < 0.5) continue;
1558 pixPtr = new AliMUONPad (mlem->GetXaxis()->GetBinCenter(j),
1559 mlem->GetYaxis()->GetBinCenter(i), 0, 0, cont1);
1560 fPixArray->Add(pixPtr);
1561 }
1562 }
c0a16418 1563}
1564
1565//_____________________________________________________________________________
1566AliMUONClusterFinderMLEM&
1567AliMUONClusterFinderMLEM::operator=(const AliMUONClusterFinderMLEM& rhs)
1568{
1569/// Protected assignement operator
1570
1571 if (this == &rhs) return *this;
1572
1573 AliFatal("Not implemented.");
1574
1575 return *this;
1576}
1577
c0a16418 1578//_____________________________________________________________________________
1579void AliMUONClusterFinderMLEM::AddVirtualPad(AliMUONCluster& cluster)
1580{
2abdae6e 1581 /// Add virtual pad (with small charge) to improve fit for clusters
1582 /// with number of pads == 2 per direction
c0a16418 1583
2abdae6e 1584 // Find out non-bending and bending planes
1585 Int_t nonb[2] = {1, 0}; // non-bending and bending cathodes
c0a16418 1586
2abdae6e 1587 TVector2 dim0 = cluster.MinPadDimensions(0, 0, kTRUE);
1588 TVector2 dim1 = cluster.MinPadDimensions(1, 0, kTRUE);
1589 if (dim0.X() < dim1.X() - fgkDistancePrecision) {
1590 nonb[0] = 0;
1591 nonb[1] = 1;
1592 }
1593
1594 Bool_t same = kFALSE;
1595 if (TMath::Abs(dim0.Y()-dim1.Y()) < fgkDistancePrecision) same = kTRUE; // the same pad size on both planes
1596
1597 AliMpIntPair cn;
1598 Bool_t check[2] = {kFALSE, kFALSE};
1599 Int_t nxy[2];
1600 nxy[0] = nxy[1] = 0;
1601 for (Int_t inb = 0; inb < 2; ++inb) {
1602 cn = cluster.NofPads(nonb[inb], 0, kTRUE);
1603 if (inb == 0 && cn.GetFirst() == 2) check[inb] = kTRUE; // check non-bending plane
1604 else if (inb == 1 && cn.GetSecond() == 2) check[inb] = kTRUE; // check bending plane
1605 if (same) {
1606 nxy[0] = TMath::Max (nxy[0], cn.GetFirst());
1607 nxy[1] = TMath::Max (nxy[1], cn.GetSecond());
1608 if (inb == 0 && nxy[0] < 2) nonb[inb] = !nonb[inb];
1609 else if (inb == 1 && cn.GetSecond() < 2) nonb[inb] = !nonb[inb];
c0a16418 1610 }
1611 }
2abdae6e 1612 if (same) {
1613 if (nxy[0] > 2) check[0] = kFALSE;
1614 if (nxy[1] > 2) check[1] = kFALSE;
1615 }
1616 if (!check[0] && !check[1]) return;
1617
1618 for (Int_t inb = 0; inb < 2; ++inb) {
1619 if (!check[inb]) continue;
1620
1621 // Find pads with maximum and next to maximum charges
1622 Int_t maxPads[2] = {-1, -1};
1623 Double_t amax[2] = {0};
1624 Int_t mult = cluster.Multiplicity();
1625 for (Int_t j = 0; j < mult; ++j) {
1626 AliMUONPad *pad = cluster.Pad(j);
1627 if (pad->Cathode() != nonb[inb]) continue;
1628 for (Int_t j2 = 0; j2 < 2; ++j2) {
1629 if (pad->Charge() > amax[j2]) {
1630 if (j2 == 0) { amax[1] = amax[0]; maxPads[1] = maxPads[0]; }
1631 amax[j2] = pad->Charge();
1632 maxPads[j2] = j;
1633 break;
1634 }
c0a16418 1635 }
1636 }
c0a16418 1637
2abdae6e 1638 // Find min and max dimensions of the cluster
1639 Double_t limits[2] = {9999, -9999};
1640 for (Int_t j = 0; j < mult; ++j) {
1641 AliMUONPad *pad = cluster.Pad(j);
1642 if (pad->Cathode() != nonb[inb]) continue;
1643 if (pad->Coord(inb) < limits[0]) limits[0] = pad->Coord(inb);
1644 if (pad->Coord(inb) > limits[1]) limits[1] = pad->Coord(inb);
c0a16418 1645 }
2abdae6e 1646
1647 // Loop over max and next to max pads
1648 Bool_t add = kFALSE;
1649 Int_t idirAdd = 0;
1650 for (Int_t j = 0; j < 2; ++j) {
1651 if (j == 1) {
1652 if (maxPads[j] < 0) continue;
1653 if (!add) break;
1654 if (amax[1] / amax[0] < 0.5) break;
c0a16418 1655 }
2abdae6e 1656 // Check if pad at the cluster limit
1657 AliMUONPad *pad = cluster.Pad(maxPads[j]);
1658 Int_t idir = 0;
1659 if (TMath::Abs(pad->Coord(inb)-limits[0]) < fgkDistancePrecision) idir = -1;
1660 else if (TMath::Abs(pad->Coord(inb)-limits[1]) < fgkDistancePrecision) idir = 1;
1661 else {
1662 //cout << " *** Pad not at the cluster limit: " << j << endl;
1663 break;
c0a16418 1664 }
2abdae6e 1665 if (j == 1 && idir == idirAdd) break; // this pad has been already added
1666
1667 // Add pad (if it exists)
1668 TVector2 pos;
1669 if (inb == 0) pos.Set (pad->X() + idir * (pad->DX()+fgkDistancePrecision), pad->Y());
1670 else pos.Set (pad->X(), pad->Y() + idir * (pad->DY()+fgkDistancePrecision));
1671 //AliMpPad mppad = fSegmentation[nonb[inb]]->PadByPosition(pos,kTRUE);
1672 AliMpPad mppad = fSegmentation[nonb[inb]]->PadByPosition(pos,kFALSE);
1673 if (!mppad.IsValid()) continue; // non-existing pad
1674 cn = mppad.GetIndices();
1675 AliMUONPad muonPad(fDetElemId, nonb[inb], cn.GetFirst(), cn.GetSecond(),
1676 mppad.Position().X(), mppad.Position().Y(),
1677 mppad.Dimensions().X(), mppad.Dimensions().Y(), 0);
1678 if (inb == 0) muonPad.SetCharge(TMath::Min (amax[j]/100, 5.));
8c718ff8 1679 //else muonPad.SetCharge(TMath::Min (amax[j]/15, fgkZeroSuppression));
1680 else muonPad.SetCharge(TMath::Min (amax[j]/15, 6.));
2abdae6e 1681 if (muonPad.Charge() < 1.) muonPad.SetCharge(1.);
1682 muonPad.SetReal(kFALSE);
1683 if (fDebug) printf(" ***** Add virtual pad in %d direction ***** %f %f %f %3d %3d %f %f \n",
1684 inb, muonPad.Charge(), muonPad.X(), muonPad.Y(), muonPad.Ix(),
1685 muonPad.Iy(), muonPad.DX(), muonPad.DY());
1686 cluster.AddPad(muonPad); // add pad to the cluster
1687 add = kTRUE;
1688 idirAdd = idir;
1689 }
1690 }
c0a16418 1691}
1692
1693//_____________________________________________________________________________
1694void AliMUONClusterFinderMLEM::PadsInXandY(AliMUONCluster& cluster,
1695 Int_t &nInX, Int_t &nInY) const
1696{
1697 /// Find number of pads in X and Y-directions (excluding virtual ones and
1698 /// overflows)
1699
8c718ff8 1700 //Int_t statusToTest = 1;
1701 Int_t statusToTest = fgkUseForFit;
c0a16418 1702
8c718ff8 1703 //if ( nInX < 0 ) statusToTest = 0;
1704 if ( nInX < 0 ) statusToTest = fgkZero;
c0a16418 1705
1706 Bool_t mustMatch(kTRUE);
1707
1708 AliMpIntPair cn = cluster.NofPads(statusToTest,mustMatch);
1709
1710 nInX = cn.GetFirst();
1711 nInY = cn.GetSecond();
1712}
1713
1714//_____________________________________________________________________________
1715void AliMUONClusterFinderMLEM::RemovePixel(Int_t i)
1716{
1717 /// Remove pixel at index i
1718 AliMUONPad* pixPtr = Pixel(i);
1719 fPixArray->RemoveAt(i);
1720 delete pixPtr;
1721}
1722
1723//_____________________________________________________________________________
1724void AliMUONClusterFinderMLEM::Simple(AliMUONCluster& cluster)
1725{
1726/// Process simple cluster (small number of pads) without EM-procedure
1727
1728 Int_t nForFit = 1, clustFit[1] = {0}, nfit;
1729 Double_t parOk[3] = {0.};
1730 TObjArray *clusters[1];
1731 clusters[0] = fPixArray;
1732
1733 AliDebug(1,Form("nPix=%d",fPixArray->GetLast()+1));
1734
2abdae6e 1735 Int_t mult = cluster.Multiplicity();
1736 for (Int_t i = 0; i < mult; ++i)
c0a16418 1737 {
1738 AliMUONPad* pad = cluster.Pad(i);
8c718ff8 1739 /*
9bbd7f60 1740 if ( pad->IsSaturated())
c0a16418 1741 {
1742 pad->SetStatus(-9);
1743 }
1744 else
1745 {
1746 pad->SetStatus(1);
1747 }
8c718ff8 1748 */
1749 if (!pad->IsSaturated()) pad->SetStatus(fgkUseForFit);
c0a16418 1750 }
1751 nfit = fSplitter->Fit(cluster,1, nForFit, clustFit, clusters, parOk, fClusterList);
1752}
1753
1754//_____________________________________________________________________________
1755AliMUONPad*
1756AliMUONClusterFinderMLEM::Pixel(Int_t i) const
1757{
1758 /// Returns pixel at index i
1759 return static_cast<AliMUONPad*>(fPixArray->UncheckedAt(i));
1760}
1761
1762//_____________________________________________________________________________
1763void
1764AliMUONClusterFinderMLEM::Print(Option_t* what) const
1765{
1766 /// printout
1767 TString swhat(what);
1768 swhat.ToLower();
1769 if ( swhat.Contains("precluster") )
1770 {
1771 if ( fPreCluster) fPreCluster->Print();
1772 }
1773}
1774
1775