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