]> git.uio.no Git - u/mrichter/AliRoot.git/blame - MUON/AliMUONClusterFinderPeakFit.cxx
In macros:
[u/mrichter/AliRoot.git] / MUON / AliMUONClusterFinderPeakFit.cxx
CommitLineData
3d4412e3 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
18//-----------------------------------------------------------------------------
19/// \class AliMUONClusterFinderPeakFit
20///
21/// Clusterizer class based on simple peak finder
22///
23/// Pre-clustering is handled by AliMUONPreClusterFinder
24/// From a precluster a pixel array is built, and its local maxima are used
25/// to get pads and make the fit with up to 3 hit candidates or compute pad
26/// centers of gravity for larger number of peaks.
27///
28/// \author Laurent Aphecetche (for the "new" C++ structure) and
29/// Alexander Zinchenko, JINR Dubna, for the hardcore of it ;-)
30//-----------------------------------------------------------------------------
31
32#include "AliMUONClusterFinderPeakFit.h"
33#include "AliMUONCluster.h"
34#include "AliMUONConstants.h"
35#include "AliMUONPad.h"
36#include "AliMUONMathieson.h"
37
38#include "AliMpDEManager.h"
39#include "AliMpPad.h"
3d4412e3 40#include "AliMpVSegmentation.h"
168e9c4d 41#include "AliMpEncodePair.h"
3d4412e3 42
43#include "AliLog.h"
44#include "AliRunLoader.h"
45//#include "AliCodeTimer.h"
46
47#include <Riostream.h>
48#include <TH2.h>
49#include <TVirtualFitter.h>
50#include <TMath.h>
51//#include <TCanvas.h>
52
53/// \cond CLASSIMP
54ClassImp(AliMUONClusterFinderPeakFit)
55/// \endcond
56
57const Double_t AliMUONClusterFinderPeakFit::fgkZeroSuppression = 6; // average zero suppression value
58//const Double_t AliMUONClusterFinderMLEM::fgkDistancePrecision = 1e-6; // (cm) used to check overlaps and so on
59const Double_t AliMUONClusterFinderPeakFit::fgkDistancePrecision = 1e-3; // (cm) used to check overlaps and so on
60const TVector2 AliMUONClusterFinderPeakFit::fgkIncreaseSize(-AliMUONClusterFinderPeakFit::fgkDistancePrecision,-AliMUONClusterFinderPeakFit::fgkDistancePrecision);
61const TVector2 AliMUONClusterFinderPeakFit::fgkDecreaseSize(AliMUONClusterFinderPeakFit::fgkDistancePrecision,AliMUONClusterFinderPeakFit::fgkDistancePrecision);
62
63// Status flags for pads
64const Int_t AliMUONClusterFinderPeakFit::fgkZero = 0x0; ///< pad "basic" state
65const Int_t AliMUONClusterFinderPeakFit::fgkMustKeep = 0x1; ///< do not kill (for pixels)
66const Int_t AliMUONClusterFinderPeakFit::fgkUseForFit = 0x10; ///< should be used for fit
67const Int_t AliMUONClusterFinderPeakFit::fgkOver = 0x100; ///< processing is over
68const Int_t AliMUONClusterFinderPeakFit::fgkModified = 0x1000; ///< modified pad charge
69const Int_t AliMUONClusterFinderPeakFit::fgkCoupled = 0x10000; ///< coupled pad
70
71namespace
72{
73 //_____________________________________________________________________________
74 Double_t Param2Coef(Int_t icand, Double_t coef, Double_t *par, Int_t nHits)
75 {
76 /// Extract hit contribution scale factor from fit parameters
77
78 //Int_t nHits = TMath::Nint(par[8]);
79 if (nHits == 1) return 1.;
80 if (nHits == 2) return icand==0 ? par[2] : TMath::Max(1.-par[2],0.);
81 if (icand == 0) return par[2];
82 if (icand == 1) return TMath::Max((1.-par[2])*par[5], 0.);
83 return TMath::Max(1.-par[2]-coef,0.);
84 }
85
86 //___________________________________________________________________________
87 void
88 FitFunction(Int_t& /*notused*/, Double_t* /*notused*/,
89 Double_t& f, Double_t* par,
90 Int_t /*notused*/)
91 {
92 /// Chi2 Function to minimize: Mathieson charge distribution in 2 dimensions
93
94 TObjArray* userObjects = static_cast<TObjArray*>(TVirtualFitter::GetFitter()->GetObjectFit());
95
96 AliMUONCluster* cluster = static_cast<AliMUONCluster*>(userObjects->At(0));
97 AliMUONMathieson* mathieson = static_cast<AliMUONMathieson*>(userObjects->At(1));
98 AliMUONClusterFinderPeakFit* finder =
99 static_cast<AliMUONClusterFinderPeakFit*>(userObjects->At(2));
100
101 f = 0.0;
102 Int_t nHits = finder->GetNMax(), npads = cluster->Multiplicity();
103 Double_t qTot = cluster->Charge(), coef = 0;
104 //if (cluster->Multiplicity(0) == 0 || cluster->Multiplicity(1) == 0) qTot *= 2.;
105
106 for ( Int_t i = 0 ; i < npads; ++i )
107 {
108 AliMUONPad* pad = cluster->Pad(i);
109 // skip pads w/ saturation or other problem(s)
110 //if ( pad->Status() ) continue;
111 if ( pad->IsSaturated() ) continue;
112 Double_t charge = 0.;
113 for (Int_t j = 0; j < nHits; ++j) {
114 // Sum over hits
115 Int_t indx = 3 * j;
116 TVector2 lowerLeft = TVector2(par[indx],par[indx+1]) - pad->Position() - pad->Dimensions();
117 TVector2 upperRight(lowerLeft + pad->Dimensions()*2.0);
118 Double_t estimatedCharge = mathieson->IntXY(lowerLeft.X(),lowerLeft.Y(),
119 upperRight.X(),upperRight.Y());
120 coef = Param2Coef(j, coef, par, nHits);
121 charge += estimatedCharge * coef;
122 }
123 charge *= qTot;
124
125 Double_t delta = charge - pad->Charge();
126 delta *= delta;
127 delta /= pad->Charge();
128 f += delta;
129 }
130 f /= (qTot/npads);
131 //cout << qTot << " " << par[0] << " " << par[1] << " " << f << endl;
132 }
133}
134
135//_____________________________________________________________________________
136AliMUONClusterFinderPeakFit::AliMUONClusterFinderPeakFit(Bool_t plot, AliMUONVClusterFinder* clusterFinder)
137 : AliMUONVClusterFinder(),
138fPreClusterFinder(clusterFinder),
139fPreCluster(0x0),
140fClusterList(),
141fMathieson(0x0),
142fEventNumber(0),
143fDetElemId(-1),
144fClusterNumber(0),
145fNMax(0),
146fHistAnode(0x0),
147fPixArray(new TObjArray(20)),
148fDebug(0),
149fPlot(plot),
150fNClusters(0),
151fNAddVirtualPads(0)
152{
153 /// Constructor
154
72dae9ff 155 fkSegmentation[1] = fkSegmentation[0] = 0x0;
3d4412e3 156
157 if (fPlot) fDebug = 1;
158}
159
160//_____________________________________________________________________________
161AliMUONClusterFinderPeakFit::~AliMUONClusterFinderPeakFit()
162{
163/// Destructor
164 delete fPixArray; fPixArray = 0;
165 delete fPreClusterFinder;
166 AliInfo(Form("Total clusters %d AddVirtualPad needed %d",
167 fNClusters,fNAddVirtualPads));
168 delete fMathieson;
169}
170
171//_____________________________________________________________________________
172Bool_t
9e41a340 173AliMUONClusterFinderPeakFit::Prepare(Int_t detElemId, TObjArray* pads[2],
110edb51 174 const AliMpArea& area,
175 const AliMpVSegmentation* seg[2])
3d4412e3 176{
177 /// Prepare for clustering
99c136e1 178// AliCodeTimerAuto("",0)
3d4412e3 179
180 for ( Int_t i = 0; i < 2; ++i )
181 {
72dae9ff 182 fkSegmentation[i] = seg[i];
3d4412e3 183 }
184
185 // Find out the DetElemId
186 fDetElemId = detElemId;
187
188 // find out current event number, and reset the cluster number
33c3c91a 189 AliRunLoader *runLoader = AliRunLoader::Instance();
7deb8eb0 190 fEventNumber = runLoader ? runLoader->GetEventNumber() : 0;
3d4412e3 191 fClusterNumber = -1;
192 fClusterList.Delete();
193
194 AliDebug(3,Form("EVT %d DE %d",fEventNumber,fDetElemId));
195
4e51cfd2 196 AliMq::Station12Type stationType = AliMpDEManager::GetStation12Type(fDetElemId);
3d4412e3 197
198 Float_t kx3 = AliMUONConstants::SqrtKx3();
199 Float_t ky3 = AliMUONConstants::SqrtKy3();
200 Float_t pitch = AliMUONConstants::Pitch();
201
4e51cfd2 202 if ( stationType == AliMq::kStation1 )
3d4412e3 203 {
204 kx3 = AliMUONConstants::SqrtKx3St1();
205 ky3 = AliMUONConstants::SqrtKy3St1();
206 pitch = AliMUONConstants::PitchSt1();
207 }
208
209 delete fMathieson;
210 fMathieson = new AliMUONMathieson;
211
212 fMathieson->SetPitch(pitch);
213 fMathieson->SetSqrtKx3AndDeriveKx2Kx4(kx3);
214 fMathieson->SetSqrtKy3AndDeriveKy2Ky4(ky3);
215
216 if ( fPreClusterFinder->NeedSegmentation() )
217 {
218 return fPreClusterFinder->Prepare(detElemId,pads,area,seg);
219 }
220 else
221 {
222 return fPreClusterFinder->Prepare(detElemId,pads,area);
223 }
224}
225
226//_____________________________________________________________________________
227AliMUONCluster*
228AliMUONClusterFinderPeakFit::NextCluster()
229{
230 /// Return next cluster
99c136e1 231// AliCodeTimerAuto("",0)
3d4412e3 232
233 // if the list of clusters is not void, pick one from there
234 TObject* o = fClusterList.At(++fClusterNumber);
235 if ( o != 0x0 ) return static_cast<AliMUONCluster*>(o);
236
237 //FIXME : at this point, must check whether we've used all the digits
238 //from precluster : if not, let the preclustering know about those unused
239 //digits, so it can reuse them
240
241 // if the cluster list is exhausted, we need to go to the next
242 // pre-cluster and treat it
243
57f3728e 244 fClusterList.Delete(); // reset the list of clusters for this pre-cluster
245 fClusterNumber = -1;
246
3d4412e3 247 fPreCluster = fPreClusterFinder->NextCluster();
248
249 if (!fPreCluster)
250 {
251 // we are done
252 return 0x0;
253 }
254
3d4412e3 255 WorkOnPreCluster();
256
257 // WorkOnPreCluster may have used only part of the pads, so we check that
258 // now, and let the unused pads be reused by the preclustering...
259
260 Int_t mult = fPreCluster->Multiplicity();
261 for ( Int_t i = 0; i < mult; ++i )
262 {
263 AliMUONPad* pad = fPreCluster->Pad(i);
264 if ( !pad->IsUsed() )
265 {
266 fPreClusterFinder->UsePad(*pad);
267 }
268 }
269
270 return NextCluster();
271}
272
273//_____________________________________________________________________________
274Bool_t
275AliMUONClusterFinderPeakFit::WorkOnPreCluster()
276{
277 /// Starting from a precluster, builds a pixel array, and then
278 /// extract clusters from this array
279
99c136e1 280 // AliCodeTimerAuto("",0)
3d4412e3 281
282 if (fDebug) {
7deb8eb0 283 cout << " *** Event # " << fEventNumber
3d4412e3 284 << " det. elem.: " << fDetElemId << endl;
285 for (Int_t j = 0; j < fPreCluster->Multiplicity(); ++j) {
286 AliMUONPad* pad = fPreCluster->Pad(j);
287 printf(" bbb %3d %1d %8.4f %8.4f %8.4f %8.4f %6.1f %3d %3d %2d %1d %1d \n",
288 j, pad->Cathode(), pad->Coord(0), pad->Coord(1), pad->DX()*2, pad->DY()*2,
289 pad->Charge(), pad->Ix(), pad->Iy(), pad->Status(), pad->IsReal(), pad->IsSaturated());
290 }
291 }
292
293 AliMUONCluster* cluster = CheckPrecluster(*fPreCluster);
294 if (!cluster) return kFALSE;
295
296 BuildPixArray(*cluster);
297
298 if ( fPixArray->GetLast() < 0 )
299 {
300 AliDebug(1,"No pixel for the above cluster");
301 delete cluster;
302 return kFALSE;
303 }
304
305 Int_t nMax = 1, localMax[100], maxPos[100] = {0};
306 Double_t maxVal[100];
307
308 nMax = FindLocalMaxima(fPixArray, localMax, maxVal); // find local maxima
309
310 if (nMax > 1) TMath::Sort(nMax, maxVal, maxPos, kTRUE); // in descending order
311
312 if (nMax <= 3) {
313 FindClusterFit(*cluster, localMax, maxPos, nMax);
314 } else {
315 for (Int_t i = 0; i < nMax; ++i)
316 {
317 FindClusterCOG(*cluster, localMax, maxPos[i]);
318 }
319 }
320
321 delete cluster;
322 if (fPlot == 0) {
323 delete fHistAnode;
324 fHistAnode = 0x0;
325 }
326 return kTRUE;
327}
328
329//_____________________________________________________________________________
330Bool_t
331AliMUONClusterFinderPeakFit::Overlap(const AliMUONPad& pad, const AliMUONPad& pixel)
332{
333 /// Check if the pad and the pixel overlaps
334
335 // make a fake pad from the pixel
336 AliMUONPad tmp(pad.DetElemId(),pad.Cathode(),pad.Ix(),pad.Iy(),
337 pixel.Coord(0),pixel.Coord(1),
338 pixel.Size(0),pixel.Size(1),0);
339
340 return AliMUONPad::AreOverlapping(pad,tmp,fgkDecreaseSize);
341}
342
343//_____________________________________________________________________________
344AliMUONCluster*
345AliMUONClusterFinderPeakFit::CheckPrecluster(const AliMUONCluster& origCluster)
346{
347 /// Check precluster in order to attempt to simplify it (mostly for
348 /// two-cathode preclusters)
349
99c136e1 350 // AliCodeTimerAuto("",0)
3d4412e3 351
352 // Disregard small clusters (leftovers from splitting or noise)
353 if ((origCluster.Multiplicity()==1 || origCluster.Multiplicity()==2) &&
110edb51 354 origCluster.Charge(0)+origCluster.Charge(1) < 1.525) // JC: adc -> fc
3d4412e3 355 {
356 return 0x0;
357 }
358
359 AliMUONCluster* cluster = new AliMUONCluster(origCluster);
360
361 AliDebug(2,"Start of CheckPreCluster=");
362 //StdoutToAliDebug(2,cluster->Print("full"));
363
364 AliMUONCluster* rv(0x0);
365
366 if (cluster->Multiplicity(0) && cluster->Multiplicity(1))
367 {
368 rv = CheckPreclusterTwoCathodes(cluster);
369 }
370 else
371 {
372 rv = cluster;
373 }
374 return rv;
375}
376
377//_____________________________________________________________________________
378AliMUONCluster*
379AliMUONClusterFinderPeakFit::CheckPreclusterTwoCathodes(AliMUONCluster* cluster)
380{
381 /// Check two-cathode cluster
382
383 Int_t npad = cluster->Multiplicity();
384 Int_t* flags = new Int_t[npad];
385 for (Int_t j = 0; j < npad; ++j) flags[j] = 0;
386
387 // Check pad overlaps
388 for ( Int_t i = 0; i < npad; ++i)
389 {
390 AliMUONPad* padi = cluster->Pad(i);
391 if ( padi->Cathode() != 0 ) continue;
392 for (Int_t j = i+1; j < npad; ++j)
393 {
394 AliMUONPad* padj = cluster->Pad(j);
395 if ( padj->Cathode() != 1 ) continue;
396 if ( !AliMUONPad::AreOverlapping(*padi,*padj,fgkDecreaseSize) ) continue;
397 flags[i] = flags[j] = 1; // mark overlapped pads
398 }
399 }
400
401 // Check if all pads overlap
402 Int_t nFlags=0;
403 for (Int_t i = 0; i < npad; ++i)
404 {
405 if (!flags[i]) ++nFlags;
406 }
407
408 if (nFlags > 0)
409 {
410 // not all pads overlap.
411 if (fDebug) cout << " nFlags: " << nFlags << endl;
412 TObjArray toBeRemoved;
413 for (Int_t i = 0; i < npad; ++i)
414 {
415 AliMUONPad* pad = cluster->Pad(i);
416 if (flags[i]) continue;
417 Int_t cath = pad->Cathode();
418 Int_t cath1 = TMath::Even(cath);
419 // Check for edge effect (missing pads on the _other_ cathode)
6e97fbb8 420 AliMpPad mpPad
421 = fkSegmentation[cath1]->PadByPosition(pad->Position().X(),pad->Position().Y(),kFALSE);
3d4412e3 422 if (!mpPad.IsValid()) continue;
423 //if (nFlags == 1 && pad->Charge() < fgkZeroSuppression * 3) continue;
110edb51 424 if (nFlags == 1 && pad->Charge() < 3.05) continue; // JC: adc -> fc
3d4412e3 425 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
426 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),pad->Charge()));
427 toBeRemoved.AddLast(pad);
428 fPreCluster->Pad(i)->Release();
429 }
430 Int_t nRemove = toBeRemoved.GetEntriesFast();
431 for ( Int_t i = 0; i < nRemove; ++i )
432 {
433 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
434 }
435 }
436
437 // Check correlations of cathode charges
438 if ( !cluster->IsSaturated() && cluster->ChargeAsymmetry() > 1 )
439 {
440 // big difference
441 Int_t cathode = cluster->MaxRawChargeCathode();
442 Int_t imin(-1);
443 Int_t imax(-1);
444 Double_t cmax(0);
445 Double_t cmin(1E9);
446
447 // get min and max pad charges on the cathode opposite to the
448 // max pad (given by MaxRawChargeCathode())
449 //
450 Int_t mult = cluster->Multiplicity();
451 for ( Int_t i = 0; i < mult; ++i )
452 {
453 AliMUONPad* pad = cluster->Pad(i);
454 if ( pad->Cathode() != cathode || !pad->IsReal() )
455 {
456 // only consider pads in the opposite cathode, and
457 // only consider real pads (i.e. exclude the virtual ones)
458 continue;
459 }
460 if ( pad->Charge() < cmin )
461 {
462 cmin = pad->Charge();
463 imin = i;
464 if (imax < 0) {
465 imax = imin;
466 cmax = cmin;
467 }
468 }
469 else if ( pad->Charge() > cmax )
470 {
471 cmax = pad->Charge();
472 imax = i;
473 }
474 }
475 AliDebug(2,Form("Pad imin,imax %d,%d cmin,cmax %e,%e",
476 imin,imax,cmin,cmax));
477 //
478 // arrange pads according to their distance to the max, normalized
479 // to the pad size
480 Double_t* dist = new Double_t[mult];
481 Double_t dxMin(1E9);
482 Double_t dyMin(1E9);
483 Double_t dmin(0);
484
485 AliMUONPad* padmax = cluster->Pad(imax);
486
487 for ( Int_t i = 0; i < mult; ++i )
488 {
489 dist[i] = 0.0;
490 if ( i == imax) continue;
491 AliMUONPad* pad = cluster->Pad(i);
492 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
493 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
494 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
495 dist[i] = TMath::Sqrt(dx*dx+dy*dy);
496 if ( i == imin )
497 {
498 dmin = dist[i] + 1E-3; // distance to the pad with minimum charge
499 dxMin = dx;
500 dyMin = dy;
501 }
502 }
503
504 TMath::Sort(mult,dist,flags,kFALSE); // in ascending order
505 Double_t xmax(-1), distPrev(999);
506 TObjArray toBeRemoved;
507
508 for ( Int_t i = 0; i < mult; ++i )
509 {
510 Int_t indx = flags[i];
511 AliMUONPad* pad = cluster->Pad(indx);
512 if ( pad->Cathode() != cathode || !pad->IsReal() ) continue;
513 if ( dist[indx] > dmin )
514 {
515 // farther than the minimum pad
516 Double_t dx = (pad->X()-padmax->X())/padmax->DX()/2.0;
517 Double_t dy = (pad->Y()-padmax->Y())/padmax->DY()/2.0;
518 dx *= dxMin;
519 dy *= dyMin;
520 if (dx >= 0 && dy >= 0) continue;
521 if (TMath::Abs(dx) > TMath::Abs(dy) && dx >= 0) continue;
522 if (TMath::Abs(dy) > TMath::Abs(dx) && dy >= 0) continue;
523 }
524 if (dist[indx] > distPrev + 1) break; // overstepping empty pads
525 if ( pad->Charge() <= cmax || TMath::Abs(dist[indx]-xmax) < 1E-3 )
526 {
527 // release pad
528 if (TMath::Abs(dist[indx]-xmax) < 1.e-3)
529 {
530 cmax = TMath::Max(pad->Charge(),cmax);
531 }
532 else
533 {
534 cmax = pad->Charge();
535 }
536 xmax = dist[indx];
537 distPrev = dist[indx];
538 AliDebug(2,Form("Releasing the following pad : de,cath,ix,iy %d,%d,%d,%d charge %e",
539 fDetElemId,pad->Cathode(),pad->Ix(),pad->Iy(),
540 pad->Charge()));
541
542 toBeRemoved.AddLast(pad);
543 fPreCluster->Pad(indx)->Release();
544 }
545 }
546 Int_t nRemove = toBeRemoved.GetEntriesFast();
547 for ( Int_t i = 0; i < nRemove; ++i )
548 {
549 cluster->RemovePad(static_cast<AliMUONPad*>(toBeRemoved.UncheckedAt(i)));
550 }
551 delete[] dist;
552 } // if ( !cluster->IsSaturated() &&
553
554 delete[] flags;
555
556 AliDebug(2,"End of CheckPreClusterTwoCathodes=");
557 //StdoutToAliDebug(2,cluster->Print("full"));
558
559 return cluster;
560}
561
562//_____________________________________________________________________________
563void
564AliMUONClusterFinderPeakFit::CheckOverlaps()
565{
566 /// For debug only : check if some pixels overlap...
567
568 Int_t nPix = fPixArray->GetLast()+1;
569 Int_t dummy(0);
570
571 for ( Int_t i = 0; i < nPix; ++i )
572 {
573 AliMUONPad* pixelI = Pixel(i);
574 AliMUONPad pi(dummy,dummy,dummy,dummy,
575 pixelI->Coord(0),pixelI->Coord(1),
576 pixelI->Size(0),pixelI->Size(1),0.0);
577
578 for ( Int_t j = i+1; j < nPix; ++j )
579 {
580 AliMUONPad* pixelJ = Pixel(j);
581 AliMUONPad pj(dummy,dummy,dummy,dummy,
582 pixelJ->Coord(0),pixelJ->Coord(1),
583 pixelJ->Size(0),pixelJ->Size(1),0.0);
584 AliMpArea area;
585
586 if ( AliMUONPad::AreOverlapping(pi,pj,fgkDecreaseSize,area) )
587 {
588 AliInfo(Form("The following 2 pixels (%d and %d) overlap !",i,j));
589 /*
590 StdoutToAliInfo(pixelI->Print();
591 cout << " Surface = " << pixelI->Size(0)*pixelI->Size(1)*4 << endl;
592 pixelJ->Print();
593 cout << " Surface = " << pixelJ->Size(0)*pixelJ->Size(1)*4 << endl;
594 cout << " Area surface = " << area.Dimensions().X()*area.Dimensions().Y()*4 << endl;
595 cout << "-------" << endl;
596 );
597 */
598 }
599 }
600 }
601}
602
603//_____________________________________________________________________________
604void AliMUONClusterFinderPeakFit::BuildPixArray(AliMUONCluster& cluster)
605{
606 /// Build pixel array
607
608 Int_t npad = cluster.Multiplicity();
609 if (npad<=0)
610 {
611 AliWarning("Got no pad at all ?!");
612 }
613
614 fPixArray->Delete();
615 BuildPixArrayOneCathode(cluster);
616
617// StdoutToAliDebug(2,cout << "End of BuildPixelArray:" << endl;
618// fPixArray->Print(););
619 //CheckOverlaps();//FIXME : this is for debug only. Remove it.
620}
621
622//_____________________________________________________________________________
623void AliMUONClusterFinderPeakFit::BuildPixArrayOneCathode(AliMUONCluster& cluster)
624{
625 /// Build the pixel array
626
627// AliDebug(2,Form("cluster.Multiplicity=%d",cluster.Multiplicity()));
628
629 TVector2 dim = cluster.MinPadDimensions (-1, kFALSE);
8c14ac2e 630 Double_t width[2] = {dim.X(), dim.Y()}, xy0[2] = { 0.0, 0.0 };
3d4412e3 631 Int_t found[2] = {0}, mult = cluster.Multiplicity();
632
633 for ( Int_t i = 0; i < mult; ++i) {
634 AliMUONPad* pad = cluster.Pad(i);
635 for (Int_t j = 0; j < 2; ++j) {
636 if (found[j] == 0 && TMath::Abs(pad->Size(j)-width[j]) < fgkDistancePrecision) {
637 xy0[j] = pad->Coord(j);
638 found[j] = 1;
639 }
640 }
641 if (found[0] && found[1]) break;
642 }
643
644 Double_t min[2], max[2];
645 Int_t cath0 = 0, cath1 = 1;
646 if (cluster.Multiplicity(0) == 0) cath0 = 1;
647 else if (cluster.Multiplicity(1) == 0) cath1 = 0;
648
6e97fbb8 649 Double_t leftDownX, leftDownY;
650 cluster.Area(cath0).LeftDownCorner(leftDownX, leftDownY);
651 Double_t rightUpX, rightUpY;
652 cluster.Area(cath0).RightUpCorner(rightUpX, rightUpY);
653 min[0] = leftDownX;
654 min[1] = leftDownY;
655 max[0] = rightUpX;
656 max[1] = rightUpY;
3d4412e3 657 if (cath1 != cath0) {
6e97fbb8 658 cluster.Area(cath1).LeftDownCorner(leftDownX, leftDownY);
659 cluster.Area(cath1).RightUpCorner(rightUpX, rightUpY);
660 min[0] = TMath::Max (min[0], leftDownX);
661 min[1] = TMath::Max (min[1], leftDownY);
662 max[0] = TMath::Min (max[0], rightUpX);
663 max[1] = TMath::Min (max[1], rightUpY);
3d4412e3 664 }
665
666 // Adjust limits
667 if (cath0 != cath1) {
668 TVector2 dim0 = cluster.MinPadDimensions (0, -1, kFALSE);
669 TVector2 dim1 = cluster.MinPadDimensions (1, -1, kFALSE);
670 if (TMath::Abs(dim0.Y()-dim1.Y()) < fgkDistancePrecision) {
671 // The same size of pads on both cathodes - check position
672 AliMUONPad* pad0 = cluster.Pad(0);
673 for ( Int_t i = 1; i < mult; ++i) {
674 AliMUONPad* pad = cluster.Pad(i);
675 if (pad->Cathode() == pad0->Cathode()) continue;
676 Double_t dist = TMath::Abs (pad0->Coord(1) - pad->Coord(1));
677 Double_t dd = dist - Int_t(dist/width[1]/2.) * width[1] * 2.;
678 if (TMath::Abs(dd/width[1]/2.-0.5) < fgkDistancePrecision) {
679 // Half pad shift between cathodes
680 width[0] /= 2.;
681 width[1] /= 2.;
682 }
683 break;
684 }
685 }
686 }
687
688 Int_t nbins[2];
689 for (Int_t i = 0; i < 2; ++i) {
690 Double_t dist = (min[i] - xy0[i]) / width[i] / 2;
691 if (TMath::Abs(dist) < 1.e-6) dist = -1.e-6;
692 min[i] = xy0[i] + (TMath::Nint(dist-TMath::Sign(1.e-6,dist))
693 + TMath::Sign(0.5,dist)) * width[i] * 2;
694 nbins[i] = TMath::Nint ((max[i] - min[i]) / width[i] / 2);
695 if (nbins[i] == 0) ++nbins[i];
696 max[i] = min[i] + nbins[i] * width[i] * 2;
697 //cout << dist << " " << min[i] << " " << max[i] << " " << nbins[i] << endl;
698 }
699
700 // Book histogram
701 TH2D *hist1 = new TH2D ("Grid", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
702 TH2D *hist2 = new TH2D ("Entries", "", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
703 TAxis *xaxis = hist1->GetXaxis();
704 TAxis *yaxis = hist1->GetYaxis();
705
706 // Fill histogram
707 for ( Int_t i = 0; i < mult; ++i) {
708 AliMUONPad* pad = cluster.Pad(i);
709 Int_t ix0 = xaxis->FindBin(pad->X());
710 Int_t iy0 = yaxis->FindBin(pad->Y());
711 PadOverHist(0, ix0, iy0, pad, hist1, hist2);
712 }
713
714 // Store pixels
715 for (Int_t i = 1; i <= nbins[0]; ++i) {
716 Double_t x = xaxis->GetBinCenter(i);
717 for (Int_t j = 1; j <= nbins[1]; ++j) {
110edb51 718 if (hist2->GetCellContent(i,j) < 0.01525) continue; // JC: adc -> fc
3d4412e3 719 if (cath0 != cath1) {
720 // Two-sided cluster
721 Double_t cont = hist2->GetCellContent(i,j);
722 if (cont < 999.) continue;
110edb51 723 if (cont-Int_t(cont/1000.)*1000. < 0.07625) continue; // JC: adc -> fc
3d4412e3 724 }
725 //if (hist2->GetCellContent(i,j) < 1.1 && cluster.Multiplicity(0) &&
726 // cluster.Multiplicity(1)) continue;
727 Double_t y = yaxis->GetBinCenter(j);
728 Double_t charge = hist1->GetCellContent(i,j);
729 AliMUONPad* pixPtr = new AliMUONPad(x, y, width[0], width[1], charge);
730 fPixArray->Add(pixPtr);
731 }
732 }
733 /*
734 if (fPixArray->GetEntriesFast() == 1) {
735 // Split pixel into 2
736 AliMUONPad* pixPtr = static_cast<AliMUONPad*> (fPixArray->UncheckedAt(0));
737 pixPtr->SetSize(0,width[0]/2.);
738 pixPtr->Shift(0,-width[0]/4.);
739 pixPtr = new AliMUONPad(pixPtr->X()+width[0], pixPtr->Y(), width[0]/2., width[1], pixPtr->Charge());
740 fPixArray->Add(pixPtr);
741 }
742 */
743 //fPixArray->Print();
744 delete hist1;
745 delete hist2;
746}
747
748//_____________________________________________________________________________
749void AliMUONClusterFinderPeakFit::PadOverHist(Int_t idir, Int_t ix0, Int_t iy0, AliMUONPad *pad,
750 TH2D *hist1, TH2D *hist2)
751{
752 /// "Span" pad over histogram in the direction idir
753
754 TAxis *axis = idir == 0 ? hist1->GetXaxis() : hist1->GetYaxis();
755 Int_t nbins = axis->GetNbins(), cath = pad->Cathode();
756 Double_t bin = axis->GetBinWidth(1), amask = TMath::Power(1000.,cath*1.);
757
758 Int_t nbinPad = (Int_t)(pad->Size(idir)/bin*2+fgkDistancePrecision) + 1; // number of bins covered by pad
759
760 for (Int_t i = 0; i < nbinPad; ++i) {
761 Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
762 if (ixy > nbins) break;
763 Double_t lowEdge = axis->GetBinLowEdge(ixy);
764 if (lowEdge + fgkDistancePrecision > pad->Coord(idir) + pad->Size(idir)) break;
765 if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
766 else {
767 // Fill histogram
768 Double_t cont = pad->Charge();
110edb51 769 if (hist2->GetCellContent(ix0, ixy) > 0.01525) // JC: adc -> fc
3d4412e3 770 cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont)
110edb51 771 + TMath::Min (TMath::Max(hist1->GetCellContent(ix0, ixy),cont)*0.1, 1.525); // JC: adc -> fc
3d4412e3 772 hist1->SetCellContent(ix0, ixy, cont);
773 hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
774 }
775 }
776
777 for (Int_t i = -1; i > -nbinPad; --i) {
778 Int_t ixy = idir == 0 ? ix0 + i : iy0 + i;
779 if (ixy < 1) break;
780 Double_t upEdge = axis->GetBinUpEdge(ixy);
781 if (upEdge - fgkDistancePrecision < pad->Coord(idir) - pad->Size(idir)) break;
782 if (idir == 0) PadOverHist(1, ixy, iy0, pad, hist1, hist2); // span in the other direction
783 else {
784 // Fill histogram
785 Double_t cont = pad->Charge();
110edb51 786 if (hist2->GetCellContent(ix0, ixy) > 0.01525) // JC: adc -> fc
3d4412e3 787 cont = TMath::Min (hist1->GetCellContent(ix0, ixy), cont)
110edb51 788 + TMath::Min (TMath::Max(hist1->GetCellContent(ix0, ixy),cont)*0.1,1.525); // JC: adc -> fc
3d4412e3 789 hist1->SetCellContent(ix0, ixy, cont);
790 hist2->SetCellContent(ix0, ixy, hist2->GetCellContent(ix0, ixy)+amask);
791 }
792 }
793}
794
795//_____________________________________________________________________________
796Int_t AliMUONClusterFinderPeakFit::FindLocalMaxima(TObjArray *pixArray, Int_t *localMax, Double_t *maxVal)
797{
798/// Find local maxima in pixel space
799
800 AliDebug(1,Form("nPix=%d",pixArray->GetLast()+1));
801
802 //TH2D *hist = NULL;
803 //delete ((TH2D*) gROOT->FindObject("anode"));
804 //if (pixArray == fPixArray) hist = (TH2D*) gROOT->FindObject("anode");
805 //else { hist = (TH2D*) gROOT->FindObject("anode1"); cout << hist << endl; }
806 //if (hist) hist->Delete();
807 delete fHistAnode;
808
809 Double_t xylim[4] = {999, 999, 999, 999};
810
811 Int_t nPix = pixArray->GetEntriesFast();
79393a89 812
813 if ( nPix <= 0 ) return 0;
814
3d4412e3 815 AliMUONPad *pixPtr = 0;
816 for (Int_t ipix = 0; ipix < nPix; ++ipix) {
817 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
818 for (Int_t i = 0; i < 4; ++i)
819 xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
820 }
821 for (Int_t i = 0; i < 4; ++i) xylim[i] -= pixPtr->Size(i/2);
822
823 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
824 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
825 if (pixArray == fPixArray) fHistAnode = new TH2D("anode","anode",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
826 else fHistAnode = new TH2D("anode1","anode1",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
827 for (Int_t ipix = 0; ipix < nPix; ++ipix) {
828 pixPtr = (AliMUONPad*) pixArray->UncheckedAt(ipix);
829 fHistAnode->Fill(pixPtr->Coord(0), pixPtr->Coord(1), pixPtr->Charge());
830 }
831// if (fDraw && pixArray == fPixArray) fDraw->DrawHist("c2", hist);
832
833 Int_t nMax = 0, indx, nxy = ny * nx;
834 Int_t *isLocalMax = new Int_t[nxy];
835 for (Int_t i = 0; i < nxy; ++i) isLocalMax[i] = 0;
836
837 for (Int_t i = 1; i <= ny; ++i) {
838 indx = (i-1) * nx;
839 for (Int_t j = 1; j <= nx; ++j) {
110edb51 840 if (fHistAnode->GetCellContent(j,i) < 0.07625) continue; // JC: adc -> fc
3d4412e3 841 //if (isLocalMax[indx+j-1] < 0) continue;
842 if (isLocalMax[indx+j-1] != 0) continue;
843 FlagLocalMax(fHistAnode, i, j, isLocalMax);
844 }
845 }
846
847 for (Int_t i = 1; i <= ny; ++i) {
848 indx = (i-1) * nx;
849 for (Int_t j = 1; j <= nx; ++j) {
850 if (isLocalMax[indx+j-1] > 0) {
851 localMax[nMax] = indx + j - 1;
852 maxVal[nMax++] = fHistAnode->GetCellContent(j,i);
85a72902 853 if (nMax > 99) break;
3d4412e3 854 }
855 }
85a72902 856 if (nMax > 99) {
857 AliError(" Too many local maxima !!!");
858 break;
859 }
3d4412e3 860 }
861 if (fDebug) cout << " Local max: " << nMax << endl;
862 delete [] isLocalMax;
863 return nMax;
864}
865
866//_____________________________________________________________________________
867void AliMUONClusterFinderPeakFit::FlagLocalMax(TH2D *hist, Int_t i, Int_t j, Int_t *isLocalMax)
868{
869/// Flag pixels (whether or not local maxima)
870
871 Int_t nx = hist->GetNbinsX();
872 Int_t ny = hist->GetNbinsY();
873 Int_t cont = TMath::Nint (hist->GetCellContent(j,i));
874 Int_t cont1 = 0, indx = (i-1)*nx+j-1, indx1 = 0, indx2 = 0;
875
876 Int_t ie = i + 2, je = j + 2;
877 for (Int_t i1 = i-1; i1 < ie; ++i1) {
878 if (i1 < 1 || i1 > ny) continue;
879 indx1 = (i1 - 1) * nx;
880 for (Int_t j1 = j-1; j1 < je; ++j1) {
881 if (j1 < 1 || j1 > nx) continue;
882 if (i == i1 && j == j1) continue;
883 indx2 = indx1 + j1 - 1;
884 cont1 = TMath::Nint (hist->GetCellContent(j1,i1));
885 if (cont < cont1) { isLocalMax[indx] = -1; return; }
886 else if (cont > cont1) isLocalMax[indx2] = -1;
887 else { // the same charge
888 isLocalMax[indx] = 1;
889 if (isLocalMax[indx2] == 0) {
890 FlagLocalMax(hist, i1, j1, isLocalMax);
891 if (isLocalMax[indx2] < 0) { isLocalMax[indx] = -1; return; }
892 else isLocalMax[indx2] = -1;
893 }
894 }
895 }
896 }
897 isLocalMax[indx] = 1; // local maximum
898}
899
900//_____________________________________________________________________________
72dae9ff 901void AliMUONClusterFinderPeakFit::FindClusterFit(AliMUONCluster& cluster, const Int_t *localMax,
902 const Int_t *maxPos, Int_t nMax)
3d4412e3 903{
904/// Fit pad charge distribution with nMax hit hypothesis
905
906 //if (cluster.Multiplicity(0) == 0 || cluster.Multiplicity(1) == 0) cout << cluster.Multiplicity(0) << " " << cluster.Multiplicity(1) << " " << cluster.Charge() << " " << cluster.Charge(0) << " " << cluster.Charge(1) << " " << endl;
907
908 //TH2D *hist = (TH2D*) gROOT->FindObject("anode");
909 Int_t nx = fHistAnode->GetNbinsX();
910 //Int_t ny = hist->GetNbinsY();
911 Double_t xmin = fHistAnode->GetXaxis()->GetXmin(); //- hist->GetXaxis()->GetBinWidth(1);
912 Double_t xmax = fHistAnode->GetXaxis()->GetXmax(); //+ hist->GetXaxis()->GetBinWidth(1);
913 Double_t ymin = fHistAnode->GetYaxis()->GetXmin(); //- hist->GetYaxis()->GetBinWidth(1);
914 Double_t ymax = fHistAnode->GetYaxis()->GetXmax(); //+ hist->GetYaxis()->GetBinWidth(1);
915
916 TVirtualFitter* fitter = TVirtualFitter::Fitter(0,nMax*3);
917 fitter->Clear("");
918 fitter->SetFCN(FitFunction);
919
920 Float_t stepX = 0.01; // cm
921 Float_t stepY = 0.01; // cm
922 Float_t stepQ = 0.01; //
923
924 Double_t args[10] = {-1.}; // disable printout
925
926 fitter->ExecuteCommand("SET PRINT",args,1);
927 fitter->ExecuteCommand("SET NOW",args,0); // no warnings
928
929 Int_t indx = 0;
930 fNMax = nMax;
931 for (Int_t i = 0; i < nMax; ++i) {
932 Int_t ic = localMax[maxPos[i]] / nx + 1;
933 Int_t jc = localMax[maxPos[i]] % nx + 1;
934 Double_t yc = fHistAnode->GetYaxis()->GetBinCenter(ic);
935 Double_t xc = fHistAnode->GetXaxis()->GetBinCenter(jc);
936 indx = 3 * i;
937 fitter->SetParameter(indx,"Hit X position",xc,stepX,xmin,xmax);
938 fitter->SetParameter(indx+1,"Hit Y position",yc,stepY,ymin,ymax);
939 fitter->SetParameter(indx+2,"Hit contribution",0.6,stepQ,0.,1.);
940 }
941 fitter->SetParameter(indx+2,"Hit contribution",0.,0.,0,0);
942 //fitter->SetParameter(8,"Number of hits",nMax,0.,0,0);
943
944 TObjArray userObjects;
945
946 userObjects.Add(&cluster);
947 userObjects.Add(fMathieson);
948 userObjects.Add(this);
949
950 fitter->SetObjectFit(&userObjects);
951
952 args[0] = 500.;
953 args[1] = 1.;
954 /*Int_t stat =*/ fitter->ExecuteCommand("MIGRAD",args,2);
955 //if (stat) { cout << " stat = " << stat << " " << fDetElemId << endl; /*exit(0);*/ }
956 //Int_t nvpar, nparx;
957 //Double_t amin, edm, errdef;
958 //fitter->GetStats(amin, edm, errdef, nvpar, nparx);
959 //cout << amin << endl;
960
961 Double_t qTot = cluster.Charge(), par[9] = {0.}, err[9] = {0.}, coef = 0.;
962 //par[8] = nMax;
963 for (Int_t j = 0; j < nMax; ++j) {
964 indx = 3 * j;
965 par[indx+2] = fitter->GetParameter(indx+2);
966 coef = Param2Coef(j, coef, par, nMax);
967 par[indx] = fitter->GetParameter(indx);
968 par[indx+1] = fitter->GetParameter(indx+1);
969 err[indx] = fitter->GetParError(indx);
970 err[indx+1] = fitter->GetParError(indx+1);
971
110edb51 972 if ( coef*qTot >= 2.135 ) // JC: adc -> fc
3d4412e3 973 {
974 AliMUONCluster* cluster1 = new AliMUONCluster(cluster);
975
976 cluster1->SetCharge(coef*qTot,coef*qTot);
977
978 cluster1->SetPosition(TVector2(par[indx],par[indx+1]),TVector2(err[indx],err[indx+1]));
979 cluster1->SetChi2(0.);
980
981 // FIXME: we miss some information in this cluster, as compared to
982 // the original AddRawCluster code.
983
984 AliDebug(2,Form("Adding RawCluster detElemId %4d mult %2d charge %5d (xl,yl)=(%9.6g,%9.6g)",
985 fDetElemId,cluster1->Multiplicity(),(Int_t)cluster1->Charge(),
986 cluster1->Position().X(),cluster1->Position().Y()));
987
988 fClusterList.Add(cluster1);
989 }
990 }
991}
992
993//_____________________________________________________________________________
994void AliMUONClusterFinderPeakFit::FindClusterCOG(AliMUONCluster& cluster,
72dae9ff 995 const Int_t *localMax, Int_t iMax)
3d4412e3 996{
997/// Find COG of pad charge distribution around local maximum \a iMax
998
999 //TH2D *hist = (TH2D*) gROOT->FindObject("anode");
1000 /* Just for check
1001 TCanvas* c = new TCanvas("Anode","Anode",800,600);
1002 c->cd();
1003 hist->Draw("lego1Fb"); // debug
1004 c->Update();
1005 Int_t tmp;
1006 cin >> tmp;
1007 */
1008 Int_t nx = fHistAnode->GetNbinsX();
1009 //Int_t ny = hist->GetNbinsY();
1010 Int_t ic = localMax[iMax] / nx + 1;
1011 Int_t jc = localMax[iMax] % nx + 1;
1012
1013 // Get min pad dimensions for the precluster
1014 Int_t nSides = 2;
1015 if (cluster.Multiplicity(0) == 0 || cluster.Multiplicity(1) == 0) nSides = 1;
1016 TVector2 dim0 = cluster.MinPadDimensions(0, -1, kFALSE);
1017 TVector2 dim1 = cluster.MinPadDimensions(1, -1, kFALSE);
1018 //Double_t width[2][2] = {{dim0.X(), dim0.Y()},{dim1.X(),dim1.Y()}};
1019 Int_t nonb[2] = {1, 0}; // coordinate index vs cathode
1020 if (nSides == 1 || dim0.X() < dim1.X() - fgkDistancePrecision) {
1021 nonb[0] = 0;
1022 nonb[1] = 1;
1023 }
1024 Bool_t samex = kFALSE, samey = kFALSE;
1025 if (TMath::Abs(dim0.X()-dim1.X()) < fgkDistancePrecision) samex = kTRUE; // the same X pad size on both planes
1026 if (TMath::Abs(dim0.Y()-dim1.Y()) < fgkDistancePrecision) samey = kTRUE; // the same Y pad size on both planes
1027
1028 // Drop all pixels from the array - pick up only the ones from the cluster
1029 //fPixArray->Delete();
1030
1031 Double_t wx = fHistAnode->GetXaxis()->GetBinWidth(1)/2;
1032 Double_t wy = fHistAnode->GetYaxis()->GetBinWidth(1)/2;
1033 Double_t yc = fHistAnode->GetYaxis()->GetBinCenter(ic);
1034 Double_t xc = fHistAnode->GetXaxis()->GetBinCenter(jc);
1035 Double_t cont = fHistAnode->GetCellContent(jc,ic);
57f3728e 1036 AliMUONPad pixel(xc, yc, wx, wy, cont);
1037 if (fDebug) pixel.Print("full");
3d4412e3 1038
1039 Int_t npad = cluster.Multiplicity();
1040
1041 // Pick up pads which overlap with the maximum pixel and find pads with the max signal
1042 Double_t qMax[2] = {0};
1043 AliMUONPad *matrix[2][3] = {{0x0,0x0,0x0},{0x0,0x0,0x0}};
1044 for (Int_t j = 0; j < npad; ++j)
1045 {
1046 AliMUONPad* pad = cluster.Pad(j);
57f3728e 1047 if ( Overlap(*pad,pixel) )
3d4412e3 1048 {
1049 if (fDebug) { cout << j << " "; pad->Print("full"); }
1050 if (pad->Charge() > qMax[pad->Cathode()]) {
1051 qMax[pad->Cathode()] = pad->Charge();
1052 matrix[pad->Cathode()][1] = pad;
1053 if (nSides == 1) matrix[!pad->Cathode()][1] = pad;
1054 }
1055 }
1056 }
1057 //if (nSides == 2 && (matrix[0][1] == 0x0 || matrix[1][1] == 0x0)) return; // ???
1058
1059 // Find neighbours of maxima to have 3 pads per direction (if possible)
1060 for (Int_t j = 0; j < npad; ++j)
1061 {
1062 AliMUONPad* pad = cluster.Pad(j);
1063 Int_t cath = pad->Cathode();
1064 if (pad == matrix[cath][1]) continue;
1065 Int_t nLoops = 3 - nSides;
1066
1067 for (Int_t k = 0; k < nLoops; ++k) {
1068 Int_t cath1 = cath;
1069 if (k) cath1 = !cath;
1070
1071 // Check the coordinate corresponding to the cathode (bending or non-bending case)
1072 Double_t dist = pad->Coord(nonb[cath1]) - matrix[cath][1]->Coord(nonb[cath1]);
1073 Double_t dir = TMath::Sign (1., dist);
1074 dist = TMath::Abs(dist) - pad->Size(nonb[cath1]) - matrix[cath][1]->Size(nonb[cath1]);
1075
1076 if (TMath::Abs(dist) < fgkDistancePrecision) {
1077 // Check the other coordinate
1078 dist = pad->Coord(!nonb[cath1]) - matrix[cath1][1]->Coord(!nonb[cath1]);
1079 if (TMath::Abs(dist) >
1080 TMath::Max(pad->Size(!nonb[cath1]), matrix[cath1][1]->Size(!nonb[cath1])) - fgkDistancePrecision) break;
1081 Int_t idir = TMath::Nint (dir);
1082 if (matrix[cath1][1+idir] == 0x0) matrix[cath1][1+idir] = pad;
1083 else if (pad->Charge() > matrix[cath1][1+idir]->Charge()) matrix[cath1][1+idir] = pad; // diff. segmentation
1084 //cout << pad->Coord(nonb[cath1]) << " " << pad->Coord(!nonb[cath1]) << " " << pad->Size(nonb[cath1]) << " " << pad->Size(!nonb[cath1]) << " " << pad->Charge() << endl ;
1085 break;
1086 }
1087 }
1088 }
1089
1090 Double_t coord[2] = {0.}, qAver = 0.;
1091 for (Int_t i = 0; i < 2; ++i) {
1092 Double_t q = 0.;
1093 Double_t coordQ = 0.;
1094 Int_t cath = matrix[i][1]->Cathode();
1095 if (i && nSides == 1) cath = !cath;
1096 for (Int_t j = 0; j < 3; ++j) {
1097 if (matrix[i][j] == 0x0) continue;
1098 Double_t dq = matrix[i][j]->Charge();
1099 q += dq;
1100 coordQ += dq * matrix[i][j]->Coord(nonb[cath]);
1101 //coordQ += (matrix[i][j]->Charge() * matrix[i][j]->Coord(nonb[cath]));
1102 }
1103 coord[cath] = coordQ / q;
1104 qAver = TMath::Max (qAver, q);
1105 }
1106
1107 //qAver = TMath::Sqrt(qAver);
110edb51 1108 if ( qAver >= 2.135 ) // JC: adc -> fc
3d4412e3 1109 {
1110
1111 AliMUONCluster* cluster1 = new AliMUONCluster(cluster);
1112
1113 cluster1->SetCharge(qAver,qAver);
1114 if (nonb[0] == 1)
1115 cluster1->SetPosition(TVector2(coord[1],coord[0]),TVector2(0.,0.));
1116 else
1117 cluster1->SetPosition(TVector2(coord[0],coord[1]),TVector2(0.,0.));
1118
1119 cluster1->SetChi2(0.);
1120
1121 // FIXME: we miss some information in this cluster, as compared to
1122 // the original AddRawCluster code.
1123
1124 AliDebug(2,Form("Adding RawCluster detElemId %4d mult %2d charge %5d (xl,yl)=(%9.6g,%9.6g)",
1125 fDetElemId,cluster1->Multiplicity(),(Int_t)cluster1->Charge(),
1126 cluster1->Position().X(),cluster1->Position().Y()));
1127
1128 fClusterList.Add(cluster1);
1129 }
1130}
1131
1132//_____________________________________________________________________________
1133AliMUONClusterFinderPeakFit&
1134AliMUONClusterFinderPeakFit::operator=(const AliMUONClusterFinderPeakFit& rhs)
1135{
1136/// Protected assignement operator
1137
1138 if (this == &rhs) return *this;
1139
1140 AliFatal("Not implemented.");
1141
1142 return *this;
1143}
1144
1145//_____________________________________________________________________________
1146void AliMUONClusterFinderPeakFit::PadsInXandY(AliMUONCluster& cluster,
1147 Int_t &nInX, Int_t &nInY) const
1148{
1149 /// Find number of pads in X and Y-directions (excluding virtual ones and
1150 /// overflows)
1151
1152 //Int_t statusToTest = 1;
1153 Int_t statusToTest = fgkUseForFit;
1154
1155 //if ( nInX < 0 ) statusToTest = 0;
1156 if ( nInX < 0 ) statusToTest = fgkZero;
1157
1158 Bool_t mustMatch(kTRUE);
1159
168e9c4d 1160 Long_t cn = cluster.NofPads(statusToTest,mustMatch);
3d4412e3 1161
168e9c4d 1162 nInX = AliMp::PairFirst(cn);
1163 nInY = AliMp::PairSecond(cn);
3d4412e3 1164}
1165
1166//_____________________________________________________________________________
1167void AliMUONClusterFinderPeakFit::RemovePixel(Int_t i)
1168{
1169 /// Remove pixel at index i
1170 AliMUONPad* pixPtr = Pixel(i);
1171 fPixArray->RemoveAt(i);
1172 delete pixPtr;
1173}
1174
1175//_____________________________________________________________________________
1176AliMUONPad*
1177AliMUONClusterFinderPeakFit::Pixel(Int_t i) const
1178{
1179 /// Returns pixel at index i
1180 return static_cast<AliMUONPad*>(fPixArray->UncheckedAt(i));
1181}
1182
1183//_____________________________________________________________________________
1184void
1185AliMUONClusterFinderPeakFit::Print(Option_t* what) const
1186{
1187 /// printout
1188 TString swhat(what);
1189 swhat.ToLower();
1190 if ( swhat.Contains("precluster") )
1191 {
1192 if ( fPreCluster) fPreCluster->Print();
1193 }
1194}
1195
1196