1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 // Clusterizer class developped by A. Zinchenko (Dubna)
21 #include <Riostream.h>
32 #include "AliMUONClusterFinderAZ.h"
33 #include "AliHeader.h"
36 #include "AliMUONChamber.h"
37 #include "AliMUONDigit.h"
38 #include "AliMUONHit.h"
39 #include "AliMUONChamber.h"
40 #include "AliMUONRawCluster.h"
41 #include "AliMUONClusterInput.h"
42 #include "AliMUONPixel.h"
44 #include "AliMUONLoader.h"
47 ClassImp(AliMUONClusterFinderAZ)
49 const Double_t AliMUONClusterFinderAZ::fgkCouplMin = 1.e-3; // threshold on coupling
50 AliMUONClusterFinderAZ* AliMUONClusterFinderAZ::fgClusterFinder = 0x0;
51 TMinuit* AliMUONClusterFinderAZ::fgMinuit = 0x0;
52 //FILE *lun1 = fopen("nxny.dat","w");
54 //_____________________________________________________________________________
55 AliMUONClusterFinderAZ::AliMUONClusterFinderAZ(Bool_t draw, Int_t iReco)
56 : AliMUONClusterFinderVS()
59 for (Int_t i=0; i<4; i++) {fHist[i] = 0;}
61 fSegmentation[1] = fSegmentation[0] = 0;
62 //AZ fgClusterFinder = 0x0;
64 if (!fgClusterFinder) fgClusterFinder = this;
65 if (!fgMinuit) fgMinuit = new TMinuit(8);
68 fPixArray = new TObjArray(20);
73 //_____________________________________________________________________________
74 AliMUONClusterFinderAZ::AliMUONClusterFinderAZ(const AliMUONClusterFinderAZ& rhs)
75 : AliMUONClusterFinderVS(rhs)
77 // Protected copy constructor
79 AliFatal("Not implemented.");
82 //_____________________________________________________________________________
83 AliMUONClusterFinderAZ::~AliMUONClusterFinderAZ()
86 delete fgMinuit; fgMinuit = 0; delete fPixArray; fPixArray = 0;
89 //_____________________________________________________________________________
90 void AliMUONClusterFinderAZ::FindRawClusters()
92 // To provide the same interface as in AliMUONClusterFinderVS
94 ResetRawClusters(); //AZ
95 EventLoop (gAlice->GetHeader()->GetEvent(), AliMUONClusterInput::Instance()->Chamber());
98 //_____________________________________________________________________________
99 void AliMUONClusterFinderAZ::EventLoop(Int_t nev=0, Int_t ch=0)
103 static Int_t nev0 = -1, ch0 = -1;
105 // Find requested event and chamber
106 if (nev < nev0) return;
107 else if (nev == nev0 && ch < ch0) return;
112 AliMUON *pMuon = (AliMUON*) gAlice->GetModule("MUON");
113 AliMUONChamber *iChamber = &(pMuon->Chamber(ch));
114 //fSegmentation[0] = iChamber->SegmentationModel(1);
115 //fSegmentation[1] = iChamber->SegmentationModel(2);
116 fResponse = iChamber->ResponseModel();
117 fSegmentation[0] = AliMUONClusterInput::Instance()->Segmentation2(0);
118 fSegmentation[1] = AliMUONClusterInput::Instance()->Segmentation2(1);
119 //AZ fResponse = AliMUONClusterInput::Instance()->Response();
121 Int_t ndigits[2]={9,9}, nShown[2]={0};
122 for (Int_t i=0; i<2; i++) {
123 for (Int_t j=0; j<fgkDim; j++) {fUsed[i][j]=kFALSE;}
127 if (ndigits[0] == nShown[0] && ndigits[1] == nShown[1]) {
131 return; // next chamber
133 Float_t xpad, ypad, zpad, zpad0;
134 Bool_t first = kTRUE;
135 if (fDebug) cout << " *** Event # " << nev << " chamber: " << ch << endl;
136 fnPads[0] = fnPads[1] = 0;
137 for (Int_t i=0; i<fgkDim; i++) {fPadIJ[1][i] = 0;}
139 for (Int_t iii = 0; iii<2; iii++) {
140 Int_t cath = TMath::Odd(iii);
141 ndigits[cath] = AliMUONClusterInput::Instance()->NDigits(cath); //AZ
142 if (!ndigits[0] && !ndigits[1]) return;
143 if (ndigits[cath] == 0) continue;
144 if (fDebug) cout << " ndigits: " << ndigits[cath] << " " << cath << endl;
149 Bool_t eEOC = kTRUE; // end-of-cluster
150 for (digit = 0; digit < ndigits[cath]; digit++) {
151 mdig = AliMUONClusterInput::Instance()->Digit(cath,digit);
153 // Find first unused pad
154 if (fUsed[cath][digit]) continue;
155 fSegmentation[cath]->GetPadC(fInput->DetElemId(),mdig->PadX(),mdig->PadY(),xpad,ypad,zpad0);
157 if (fUsed[cath][digit]) continue;
158 fSegmentation[cath]->GetPadC(fInput->DetElemId(),mdig->PadX(),mdig->PadY(),xpad,ypad,zpad);
159 if (TMath::Abs(zpad-zpad0) > 0.1) continue; // different slats
160 // Find a pad overlapping with the cluster
161 if (!Overlap(cath,mdig)) continue;
163 // Add pad - recursive call
165 //AZ !!!!!! Temporary fix of St1 overlap regions !!!!!!!!
166 if (cath && ch < 2) {
167 Int_t npads = fnPads[0] + fnPads[1] - 1;
168 Int_t cath1 = fPadIJ[0][npads];
169 Int_t idig = TMath::Nint (fXyq[5][npads]);
170 mdig = AliMUONClusterInput::Instance()->Digit(cath1,idig);
171 fSegmentation[cath1]->GetPadC(fInput->DetElemId(),mdig->PadX(),mdig->PadY(),xpad,ypad,zpad);
172 if (TMath::Abs(zpad-zpad0) > 0.1) zpad0 = zpad;
175 if (digit >= 0) break;
178 // No more unused pads
179 if (cath == 0) continue; // on cathode #0 - check #1
180 else return; // No more clusters
182 if (eEOC) break; // cluster found
184 if (fDebug) cout << " nPads: " << fnPads[cath] << " " << nShown[cath]+fnPads[cath] << " " << cath << endl;
185 } // for (Int_t iii = 0;
188 if (fDraw) DrawCluster(nev0, ch0);
190 // Use MLEM for cluster finder
191 Int_t nMax = 1, localMax[100], maxPos[100];
192 Double_t maxVal[100];
194 if (CheckPrecluster(nShown)) {
196 if (fnPads[0]+fnPads[1] > 50) nMax = FindLocalMaxima(localMax, maxVal);
197 if (nMax > 1) TMath::Sort(nMax, maxVal, maxPos, kTRUE); // in decreasing order
198 Int_t iSimple = 0, nInX = -1, nInY;
199 PadsInXandY(nInX, nInY);
200 if (fDebug) cout << "Pads in X and Y: " << nInX << " " << nInY << endl;
201 if (nMax == 1 && nInX < 4 && nInY < 4) iSimple = 0; //1; // simple cluster
202 for (Int_t i=0; i<nMax; i++) {
203 if (nMax > 1) FindCluster(localMax, maxPos[i]);
204 if (!MainLoop(iSimple)) cout << " MainLoop failed " << endl;
206 for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
207 if (fPadIJ[1][j] == 0) continue; // pad charge was not modified
209 fXyq[2][j] = fXyq[6][j]; // use backup charge value
214 if (fReco || Next(nev0, ch0)) goto next;
217 //_____________________________________________________________________________
218 void AliMUONClusterFinderAZ::DrawCluster(Int_t nev0, Int_t ch0)
224 Double_t p1[3]={0}, p2[3];
226 c1 = new TCanvas("c1","Clusters",0,0,600,700);
227 //c1->SetFillColor(10);
229 new TCanvas("c2","Mlem",700,0,600,350);
231 c1 = (TCanvas*) gROOT->GetListOfCanvases()->FindObject("c1");
236 // Get pointer to Alice detectors
237 AliMUON *muon = (AliMUON*) gAlice->GetModule("MUON");
241 AliRunLoader *rl = AliRunLoader::GetRunLoader();
242 AliLoader *gime = rl->GetLoader("MUONLoader");
243 AliMUONData *data = ((AliMUONLoader*)gime)->GetMUONData();
245 gime->LoadHits("READ");
246 TTree *treeH = gime->TreeH();
247 ntracks = (Int_t) treeH->GetEntries();
248 cout << " nev " << nev0 << " ntracks " << ntracks << endl;
249 gime->LoadRecPoints("READ");
250 TTree *treeR = data->TreeR();
252 data->SetTreeAddress("RC");
253 data->GetRawClusters();
259 for (Int_t cath = 0; cath<2; cath++) {
261 if (fHist[cath*2]) {fHist[cath*2]->Delete(); fHist[cath*2] = 0;}
262 if (fHist[cath*2+1]) {fHist[cath*2+1]->Delete(); fHist[cath*2+1] = 0;}
263 if (fnPads[cath] == 0) continue; // cluster on one cathode only
264 Float_t wxMin=999, wxMax=0, wyMin=999, wyMax=0;
265 Int_t minDx=0, maxDx=0, minDy=0, maxDy=0;
266 for (Int_t i=0; i<fnPads[0]+fnPads[1]; i++) {
267 if (fPadIJ[0][i] != cath) continue;
268 if (fXyq[3][i] < wxMin) {wxMin = fXyq[3][i]; minDx = i;}
269 if (fXyq[3][i] > wxMax) {wxMax = fXyq[3][i]; maxDx = i;}
270 if (fXyq[4][i] < wyMin) {wyMin = fXyq[4][i]; minDy = i;}
271 if (fXyq[4][i] > wyMax) {wyMax = fXyq[4][i]; maxDy = i;}
273 cout << minDx << maxDx << minDy << maxDy << endl;
274 Int_t nx, ny, padSize;
275 Float_t xmin=9999, xmax=-9999, ymin=9999, ymax=-9999;
276 if (TMath::Nint(fXyq[3][minDx]*1000) == TMath::Nint(fXyq[3][maxDx]*1000) &&
277 TMath::Nint(fXyq[4][minDy]*1000) == TMath::Nint(fXyq[4][maxDy]*1000)) {
278 // the same segmentation
279 cout << " Same" << endl;
280 cout << fXyq[3][minDx] << " " << fXyq[3][maxDx] << " " << fXyq[4][minDy] << " " << fXyq[4][maxDy] << endl;
281 for (Int_t i=0; i<fnPads[0]+fnPads[1]; i++) {
282 if (fPadIJ[0][i] != cath) continue;
283 if (fXyq[0][i] < xmin) xmin = fXyq[0][i];
284 if (fXyq[0][i] > xmax) xmax = fXyq[0][i];
285 if (fXyq[1][i] < ymin) ymin = fXyq[1][i];
286 if (fXyq[1][i] > ymax) ymax = fXyq[1][i];
288 xmin -= fXyq[3][minDx]; xmax += fXyq[3][minDx];
289 ymin -= fXyq[4][minDy]; ymax += fXyq[4][minDy];
290 nx = TMath::Nint ((xmax-xmin)/wxMin/2);
291 ny = TMath::Nint ((ymax-ymin)/wyMin/2);
292 cout << xmin << " " << xmax << " " << nx << " " << ymin << " " << ymax << " " << ny << endl;
293 sprintf(hName,"h%d",cath*2);
294 fHist[cath*2] = new TH2F(hName,"cluster",nx,xmin,xmax,ny,ymin,ymax);
295 //cout << fHist[cath*2] << " " << fnPads[cath] << endl;
296 for (Int_t i=0; i<fnPads[0]+fnPads[1]; i++) {
297 if (fPadIJ[0][i] != cath) continue;
298 fHist[cath*2]->Fill(fXyq[0][i],fXyq[1][i],fXyq[2][i]);
299 //cout << fXyq[0][i] << fXyq[1][i] << fXyq[2][i] << endl;
302 // different segmentation in the cluster
303 cout << " Different" << endl;
304 cout << fXyq[3][minDx] << " " << fXyq[3][maxDx] << " " << fXyq[4][minDy] << " " << fXyq[4][maxDy] << endl;
306 Int_t indx, locMin, locMax;
307 if (TMath::Nint(fXyq[3][minDx]*1000) != TMath::Nint(fXyq[3][maxDx]*1000)) {
308 // different segmentation along x
313 // different segmentation along y
319 for (Int_t i=0; i<2; i++) {
320 // loop over different pad sizes
321 if (i>0) loc = locMax;
322 padSize = TMath::Nint(fXyq[indx+3][loc]*1000);
323 xmin = 9999; xmax = -9999; ymin = 9999; ymax = -9999;
324 for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
325 if (fPadIJ[0][j] != cath) continue;
326 if (TMath::Nint(fXyq[indx+3][j]*1000) != padSize) continue;
328 xmin = TMath::Min (xmin,fXyq[0][j]);
329 xmax = TMath::Max (xmax,fXyq[0][j]);
330 ymin = TMath::Min (ymin,fXyq[1][j]);
331 ymax = TMath::Max (ymax,fXyq[1][j]);
333 xmin -= fXyq[3][loc]; xmax += fXyq[3][loc];
334 ymin -= fXyq[4][loc]; ymax += fXyq[4][loc];
335 nx = TMath::Nint ((xmax-xmin)/fXyq[3][loc]/2);
336 ny = TMath::Nint ((ymax-ymin)/fXyq[4][loc]/2);
337 sprintf(hName,"h%d",cath*2+i);
338 fHist[cath*2+i] = new TH2F(hName,"cluster",nx,xmin,xmax,ny,ymin,ymax);
339 for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
340 if (fPadIJ[0][j] != cath) continue;
341 if (TMath::Nint(fXyq[indx+3][j]*1000) != padSize) continue;
342 fHist[cath*2+i]->Fill(fXyq[0][j],fXyq[1][j],fXyq[2][j]);
345 if (nOK != fnPads[cath]) cout << " *** Too many segmentations: nPads, nOK " << fnPads[cath] << " " << nOK << endl;
346 } // if (TMath::Nint(fXyq[3][minDx]*1000)
347 } // for (Int_t cath = 0;
349 // Draw histograms and coordinates
350 for (Int_t cath=0; cath<2; cath++) {
351 if (cath == 0) ModifyHistos();
352 if (fnPads[cath] == 0) continue; // cluster on one cathode only
357 Double_t x, y, x0, y0, r1=999, r2=0;
358 if (fHist[cath*2+1]) {
360 x0 = fHist[cath*2]->GetXaxis()->GetXmin() - 1000*TMath::Cos(30*TMath::Pi()/180);
361 y0 = fHist[cath*2]->GetYaxis()->GetXmin() - 1000*TMath::Sin(30*TMath::Pi()/180);
364 for (Int_t iy=1; iy<=fHist[ihist]->GetNbinsY(); iy++) {
365 y = fHist[ihist]->GetYaxis()->GetBinCenter(iy)
366 + fHist[ihist]->GetYaxis()->GetBinWidth(iy);
367 for (Int_t ix=1; ix<=fHist[ihist]->GetNbinsX(); ix++) {
368 if (fHist[ihist]->GetCellContent(ix,iy) > 0.1) {
369 x = fHist[ihist]->GetXaxis()->GetBinCenter(ix)
370 + fHist[ihist]->GetXaxis()->GetBinWidth(ix);
371 r1 = TMath::Max (r1,TMath::Sqrt((x-x0)*(x-x0)+(y-y0)*(y-y0)));
376 for (Int_t iy=1; iy<=fHist[ihist]->GetNbinsY(); iy++) {
377 y = fHist[ihist]->GetYaxis()->GetBinCenter(iy)
378 + fHist[ihist]->GetYaxis()->GetBinWidth(iy);
379 for (Int_t ix=1; ix<=fHist[ihist]->GetNbinsX(); ix++) {
380 if (fHist[ihist]->GetCellContent(ix,iy) > 0.1) {
381 x = fHist[ihist]->GetXaxis()->GetBinCenter(ix)
382 + fHist[ihist]->GetXaxis()->GetBinWidth(ix);
383 r2 = TMath::Max (r2,TMath::Sqrt((x-x0)*(x-x0)+(y-y0)*(y-y0)));
387 cout << r1 << " " << r2 << endl;
388 } // if (fHist[cath*2+1])
390 //fHist[cath*2]->Draw("lego1");
391 fHist[cath*2]->Draw("lego1Fb");
392 //if (fHist[cath*2+1]) fHist[cath*2+1]->Draw("lego1SameAxisBb");
393 if (fHist[cath*2+1]) fHist[cath*2+1]->Draw("lego1SameAxisBbFb");
395 //fHist[cath*2+1]->Draw("lego1");
396 fHist[cath*2+1]->Draw("lego1Fb");
397 //fHist[cath*2]->Draw("lego1SameAxisBb");
398 fHist[cath*2]->Draw("lego1SameAxisFbBb");
402 } // for (Int_t cath = 0;
404 // Draw generated hits
406 hist = fHist[0] ? fHist[0] : fHist[2];
407 p2[2] = hist->GetMaximum();
409 if (c1) view = c1->Pad()->GetView();
410 cout << " *** GEANT hits *** " << endl;
413 for (Int_t i=0; i<ntracks; i++) {
415 for (AliMUONHit* mHit=(AliMUONHit*)muon->FirstHit(-1);
417 mHit=(AliMUONHit*)muon->NextHit()) {
418 if (mHit->Chamber() != ch0+1) continue; // chamber number
419 if (TMath::Abs(mHit->Z()-fZpad) > 1) continue; // different slat
420 p2[0] = p1[0] = mHit->X(); // x-pos of hit
421 p2[1] = p1[1] = mHit->Y(); // y-pos
422 if (p1[0] < hist->GetXaxis()->GetXmin() ||
423 p1[0] > hist->GetXaxis()->GetXmax()) continue;
424 if (p1[1] < hist->GetYaxis()->GetXmin() ||
425 p1[1] > hist->GetYaxis()->GetXmax()) continue;
426 // Check if track comes thru pads with signal
428 for (Int_t ihist=0; ihist<4; ihist++) {
429 if (!fHist[ihist]) continue;
430 ix = fHist[ihist]->GetXaxis()->FindBin(p1[0]);
431 iy = fHist[ihist]->GetYaxis()->FindBin(p1[1]);
432 if (fHist[ihist]->GetCellContent(ix,iy) > 0.5) {iok = 1; break;}
435 gStyle->SetLineColor(1);
436 if (TMath::Abs((Int_t)mHit->Particle()) == 13) {
437 gStyle->SetLineColor(4);
439 fxyMu[fnMu][0] = p1[0];
440 fxyMu[fnMu++][1] = p1[1];
443 if (fDebug) printf(" X=%10.4f, Y=%10.4f, Z=%10.4f\n",p1[0],p1[1],mHit->Z());
445 // Take into account track angles
446 p2[0] += mHit->Tlength() * TMath::Sin(mHit->Theta()/180*TMath::Pi())
447 * TMath::Cos(mHit->Phi()/180*TMath::Pi()) / 2;
448 p2[1] += mHit->Tlength() * TMath::Sin(mHit->Theta()/180*TMath::Pi())
449 * TMath::Sin(mHit->Phi()/180*TMath::Pi()) / 2;
450 view->WCtoNDC(p1, &xNDC[0]);
451 view->WCtoNDC(p2, &xNDC[3]);
452 for (Int_t ipad=1; ipad<3; ipad++) {
454 //c1->DrawLine(xpad[0],xpad[1],xpad[3],xpad[4]);
455 line[nLine] = new TLine(xNDC[0],xNDC[1],xNDC[3],xNDC[4]);
456 line[nLine++]->Draw();
459 } // for (AliMUONHit* mHit=
460 } // for (Int_t i=0; i<ntracks;
462 // Draw reconstructed coordinates
463 TClonesArray *listMUONrawclust = data->RawClusters(ch0);
464 AliMUONRawCluster *mRaw;
465 gStyle->SetLineColor(3);
466 cout << " *** Reconstructed hits *** " << endl;
467 if (listMUONrawclust) {
468 for (Int_t i=0; i<listMUONrawclust ->GetEntries(); i++) {
469 mRaw = (AliMUONRawCluster*)listMUONrawclust ->UncheckedAt(i);
470 if (TMath::Abs(mRaw->GetZ(0)-fZpad) > 1) continue; // different slat
471 p2[0] = p1[0] = mRaw->GetX(0); // x-pos of hit
472 p2[1] = p1[1] = mRaw->GetY(0); // y-pos
473 if (p1[0] < hist->GetXaxis()->GetXmin() ||
474 p1[0] > hist->GetXaxis()->GetXmax()) continue;
475 if (p1[1] < hist->GetYaxis()->GetXmin() ||
476 p1[1] > hist->GetYaxis()->GetXmax()) continue;
478 treeD->GetEvent(cath);
479 cout << mRaw->fMultiplicity[0] << mRaw->fMultiplicity[1] << endl;
480 for (Int_t j=0; j<mRaw->fMultiplicity[cath]; j++) {
481 Int_t digit = mRaw->fIndexMap[j][cath];
482 cout << ((AliMUONDigit*)fMuonDigits->UncheckedAt(digit))->Signal() << endl;
485 // Check if track comes thru pads with signal
487 for (Int_t ihist=0; ihist<4; ihist++) {
488 if (!fHist[ihist]) continue;
489 ix = fHist[ihist]->GetXaxis()->FindBin(p1[0]);
490 iy = fHist[ihist]->GetYaxis()->FindBin(p1[1]);
491 if (fHist[ihist]->GetCellContent(ix,iy) > 0.5) {iok = 1; break;}
494 if (fDebug) printf(" X=%10.4f, Y=%10.4f, Z=%10.4f\n",p1[0],p1[1],mRaw->GetZ(0));
496 view->WCtoNDC(p1, &xNDC[0]);
497 view->WCtoNDC(p2, &xNDC[3]);
498 for (Int_t ipad=1; ipad<3; ipad++) {
500 line[nLine] = new TLine(xNDC[0],xNDC[1],xNDC[3],xNDC[4]);
501 line[nLine++]->Draw();
504 } // for (Int_t i=0; i<listMUONrawclust ->GetEntries();
505 } // if (listMUONrawclust)
509 //_____________________________________________________________________________
510 Int_t AliMUONClusterFinderAZ::Next(Int_t &nev0, Int_t &ch0)
515 //lun = fopen("pull.dat","w");
517 for (Int_t i=0; i<fnMu; i++) {
518 // Check again if muon come thru the used pads (due to extra splitting)
519 for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
520 if (TMath::Abs(fxyMu[i][0]-fXyq[0][j])<fXyq[3][j] &&
521 TMath::Abs(fxyMu[i][1]-fXyq[1][j])<fXyq[4][j]) {
522 if (fDebug) printf("%12.3e %12.3e %12.3e %12.3e\n",fxyMu[i][2],fxyMu[i][3],fxyMu[i][4],fxyMu[i][5]);
523 if (lun) fprintf(lun,"%4d %2d %12.3e %12.3e %12.3e %12.3e\n",nev0,ch0,fxyMu[i][2],fxyMu[i][3],fxyMu[i][4],fxyMu[i][5]);
527 } // for (Int_t i=0; i<fnMu;
531 cout << " What is next? " << endl;
534 if (command[0] == 'n' || command[0] == 'N') {nev0++; ch0 = 0; } // next event
535 else if (command[0] == 'q' || command[0] == 'Q') { if (lun) fclose(lun); } // exit display
536 else if (command[0] == 'c' || command[0] == 'C') sscanf(command+1,"%d",&ch0); // new chamber
537 else if (command[0] == 'e' || command[0] == 'E') { sscanf(command+1,"%d",&nev0); ch0 = 0; } // new event
538 else return 1; // Next precluster
543 //_____________________________________________________________________________
544 void AliMUONClusterFinderAZ::ModifyHistos(void)
546 // Modify histograms to bring them to (approximately) the same size
548 Float_t hlim[4][4], hbin[4][4]; // first index - xmin, xmax, ymin, ymax
550 Float_t binMin[4] = {999,999,999,999};
552 for (Int_t i = 0; i < 4; i++) {
554 hlim[0][i] = hlim[2][i] = 999;
555 hlim[1][i] = hlim[3][i] = -999;
558 hlim[0][i] = fHist[i]->GetXaxis()->GetXmin(); // xmin
559 hlim[1][i] = fHist[i]->GetXaxis()->GetXmax(); // xmax
560 hlim[2][i] = fHist[i]->GetYaxis()->GetXmin(); // ymin
561 hlim[3][i] = fHist[i]->GetYaxis()->GetXmax(); // ymax
562 hbin[0][i] = hbin[1][i] = fHist[i]->GetXaxis()->GetBinWidth(1);
563 hbin[2][i] = hbin[3][i] = fHist[i]->GetYaxis()->GetBinWidth(1);
564 binMin[0] = TMath::Min(binMin[0],hbin[0][i]);
565 binMin[2] = TMath::Min(binMin[2],hbin[2][i]);
568 binMin[1] = binMin[0];
569 binMin[3] = binMin[2];
570 cout << " Nhist: " << nhist << endl;
572 // Adjust histo limits for cathode with different segmentation
573 for (Int_t i = 0; i < 4; i+=2) {
574 if (!fHist[i+1]) continue;
575 Int_t imin, imax, i1 = i + 1;
576 for (Int_t lim = 0; lim < 4; lim++) {
578 if (hlim[lim][i] < hlim[lim][i1]) {
585 if (TMath::Abs(hlim[lim][imin]-hlim[lim][imax])<0.01*binMin[lim]) break;
586 if (lim == 0 || lim == 2) {
588 hlim[lim][imax] -= hbin[lim][imax];
591 hlim[lim][imin] += hbin[lim][imin];
598 Int_t imnmx = 0, nExtra = 0;
599 for (Int_t lim = 0; lim < 4; lim++) {
600 if (lim == 0 || lim == 2) imnmx = TMath::LocMin(4,hlim[lim]); // find lower limit
601 else imnmx = TMath::LocMax(4,hlim[lim]); // find upper limit
603 // Adjust histogram limit
604 for (Int_t i = 0; i < 4; i++) {
605 if (!fHist[i]) continue;
606 nExtra = TMath::Nint ((hlim[lim][imnmx]-hlim[lim][i]) / hbin[lim][i]);
607 hlim[lim][i] += nExtra * hbin[lim][i];
611 // Rebuild histograms
614 Double_t x, y, cont, cmax=0;
616 for (Int_t ihist=0; ihist<4; ihist++) {
617 if (!fHist[ihist]) continue;
618 nx = TMath::Nint((hlim[1][ihist]-hlim[0][ihist])/hbin[0][ihist]);
619 ny = TMath::Nint((hlim[3][ihist]-hlim[2][ihist])/hbin[2][ihist]);
620 cout << ihist << " " << hlim[0][ihist] << " " << hlim[1][ihist] << " " << nx;
621 cout << " " << hlim[2][ihist] << " " << hlim[3][ihist] << " " << ny << endl;
622 sprintf(hName,"hh%d",ihist);
623 hist = new TH2F(hName,"hist",nx,hlim[0][ihist],hlim[1][ihist],ny,hlim[2][ihist],hlim[3][ihist]);
624 for (Int_t i=1; i<=fHist[ihist]->GetNbinsX(); i++) {
625 x = fHist[ihist]->GetXaxis()->GetBinCenter(i);
626 for (Int_t j=1; j<=fHist[ihist]->GetNbinsY(); j++) {
627 y = fHist[ihist]->GetYaxis()->GetBinCenter(j);
628 cont = fHist[ihist]->GetCellContent(i,j);
629 hist->Fill(x,y,cont);
632 cmax = TMath::Max (cmax,hist->GetMaximum());
633 sprintf(hName,"%s%d",fHist[ihist]->GetName(),ihist);
634 fHist[ihist]->Delete();
635 fHist[ihist] = new TH2F(*hist);
636 fHist[ihist]->SetName(hName);
637 fHist[ihist]->SetNdivisions(505,"Z");
640 if (fDebug) printf("%f \n",cmax);
642 for (Int_t ihist=0; ihist<4; ihist++) {
643 if (!fHist[ihist]) continue;
644 fHist[ihist]->SetMaximum(cmax);
645 fHist[ihist]->SetMinimum(0);
649 //_____________________________________________________________________________
650 void AliMUONClusterFinderAZ::AddPad(Int_t cath, Int_t digit)
652 // Add pad to the cluster
653 //AZ AliMUONDigit *mdig = (AliMUONDigit*)fMuonDigits->UncheckedAt(digit);
654 AliMUONDigit *mdig = AliMUONClusterInput::Instance()->Digit(cath,digit); //AZ
656 Int_t charge = mdig->Signal();
657 // get the center of the pad
658 Float_t xpad, ypad, zpad0; //, zpad;
659 fSegmentation[cath]->GetPadC(fInput->DetElemId(),mdig->PadX(),mdig->PadY(),xpad,ypad,zpad0);
661 Int_t isec = fSegmentation[cath]->Sector(fInput->DetElemId(),mdig->PadX(), mdig->PadY());
662 Int_t nPads = fnPads[0] + fnPads[1];
663 fXyq[0][nPads] = xpad;
664 fXyq[1][nPads] = ypad;
665 fXyq[2][nPads] = charge;
666 fXyq[3][nPads] = fSegmentation[cath]->Dpx(fInput->DetElemId(),isec)/2;
667 fXyq[4][nPads] = fSegmentation[cath]->Dpy(fInput->DetElemId(),isec)/2;
668 fXyq[5][nPads] = digit;
670 fPadIJ[0][nPads] = cath;
671 fPadIJ[1][nPads] = 0;
672 fUsed[cath][digit] = kTRUE;
673 if (fDebug) printf(" bbb %d %d %f %f %f %f %f %d\n", nPads, cath, xpad, ypad, zpad0, fXyq[3][nPads]*2, fXyq[4][nPads]*2, charge);
677 Int_t nn, ix, iy, xList[10], yList[10];
680 //AZ Int_t ndigits = fMuonDigits->GetEntriesFast();
681 Int_t ndigits = AliMUONClusterInput::Instance()->NDigits(cath); //AZ
682 fSegmentation[cath]->Neighbours(fInput->DetElemId(),mdig->PadX(),mdig->PadY(),&nn,xList,yList);
683 for (Int_t in=0; in<nn; in++) {
686 for (Int_t digit1 = 0; digit1 < ndigits; digit1++) {
687 if (digit1 == digit) continue;
688 //AZ mdig1 = (AliMUONDigit*)fMuonDigits->UncheckedAt(digit1);
689 mdig1 = AliMUONClusterInput::Instance()->Digit(cath,digit1); //AZ
690 //AZobsolete if (mdig1->Cathode() != cath) continue;
691 if (!fUsed[cath][digit1] && mdig1->PadX() == ix && mdig1->PadY() == iy) {
692 //AZ--- temporary fix on edges
693 //fSegmentation[cath]->GetPadC(mdig1->PadX(), mdig1->PadY(), xpad, ypad, zpad);
694 //if (TMath::Abs(zpad-zpad0) > 0.5) continue;
696 fUsed[cath][digit1] = kTRUE;
697 // Add pad - recursive call
700 } //for (Int_t digit1 = 0;
701 } // for (Int_t in=0;
704 //_____________________________________________________________________________
705 Bool_t AliMUONClusterFinderAZ::Overlap(Int_t cath, AliMUONDigit *mdig)
707 // Check if the pad from one cathode overlaps with a pad
708 // in the precluster on the other cathode
710 Float_t xpad, ypad, zpad;
711 fSegmentation[cath]->GetPadC(fInput->DetElemId(),mdig->PadX(),mdig->PadY(),xpad,ypad,zpad);
712 Int_t isec = fSegmentation[cath]->Sector(fInput->DetElemId(),mdig->PadX(), mdig->PadY());
714 Float_t xy1[4], xy12[4];
715 xy1[0] = xpad - fSegmentation[cath]->Dpx(fInput->DetElemId(),isec)/2;
716 xy1[1] = xy1[0] + fSegmentation[cath]->Dpx(fInput->DetElemId(),isec);
717 xy1[2] = ypad - fSegmentation[cath]->Dpy(fInput->DetElemId(),isec)/2;
718 xy1[3] = xy1[2] + fSegmentation[cath]->Dpy(fInput->DetElemId(),isec);
719 //cout << " ok " << fnPads[0]+fnPads[1] << xy1[0] << xy1[1] << xy1[2] << xy1[3] << endl;
721 Int_t cath1 = TMath::Even(cath);
722 for (Int_t i=0; i<fnPads[0]+fnPads[1]; i++) {
723 if (fPadIJ[0][i] != cath1) continue;
724 if (Overlap(xy1, i, xy12, 0)) return kTRUE;
729 //_____________________________________________________________________________
730 Bool_t AliMUONClusterFinderAZ::Overlap(Float_t *xy1, Int_t iPad, Float_t *xy12, Int_t iSkip)
732 // Check if the pads xy1 and iPad overlap and return overlap area
735 xy2[0] = fXyq[0][iPad] - fXyq[3][iPad];
736 xy2[1] = fXyq[0][iPad] + fXyq[3][iPad];
737 if (xy1[0] > xy2[1]-1.e-4 || xy1[1] < xy2[0]+1.e-4) return kFALSE;
738 xy2[2] = fXyq[1][iPad] - fXyq[4][iPad];
739 xy2[3] = fXyq[1][iPad] + fXyq[4][iPad];
740 if (xy1[2] > xy2[3]-1.e-4 || xy1[3] < xy2[2]+1.e-4) return kFALSE;
741 if (!iSkip) return kTRUE; // just check overlap (w/out computing the area)
742 xy12[0] = TMath::Max (xy1[0],xy2[0]);
743 xy12[1] = TMath::Min (xy1[1],xy2[1]);
744 xy12[2] = TMath::Max (xy1[2],xy2[2]);
745 xy12[3] = TMath::Min (xy1[3],xy2[3]);
749 //_____________________________________________________________________________
751 Bool_t AliMUONClusterFinderAZ::Overlap(Int_t i, Int_t j, Float_t *xy12, Int_t iSkip)
753 // Check if the pads i and j overlap and return overlap area
755 Float_t xy1[4], xy2[4];
756 return Overlap(xy1, xy2, xy12, iSkip);
759 //_____________________________________________________________________________
760 Bool_t AliMUONClusterFinderAZ::CheckPrecluster(Int_t *nShown)
762 // Check precluster in order to attempt to simplify it (mostly for
763 // two-cathode preclusters)
765 Int_t i1, i2, cath=0, digit=0;
766 Float_t xy1[4], xy12[4];
768 Int_t npad = fnPads[0] + fnPads[1];
770 // Disregard one-pad clusters (leftovers from splitting)
771 nShown[0] += fnPads[0];
772 nShown[1] += fnPads[1];
776 // If pads have the same size take average of pads on both cathodes
777 Int_t sameSize = (fnPads[0] && fnPads[1]) ? 1 : 0;
779 Double_t xSize = -1, ySize = 0;
780 for (Int_t i=0; i<npad; i++) {
781 if (fXyq[2][i] < 0) continue;
782 if (xSize < 0) { xSize = fXyq[3][i]; ySize = fXyq[4][i]; }
783 if (TMath::Abs(xSize-fXyq[3][i]) > 1.e-4 || TMath::Abs(ySize-fXyq[4][i]) > 1.e-4) { sameSize = 0; break; }
786 if (sameSize && fnPads[0] == 1 && fnPads[1] == 1) sameSize = 0; //AZ
787 if (sameSize && (fnPads[0] >= 2 || fnPads[1] >= 2)) {
788 nShown[0] += fnPads[0];
789 nShown[1] += fnPads[1];
790 fnPads[0] = fnPads[1] = 0;
792 for (Int_t i=0; i<npad; i++) {
793 if (fXyq[2][i] < 0) continue; // used pad
794 fXyq[2][fnPads[0]] = fXyq[2][i];
797 for (Int_t j=i+1; j<npad; j++) {
798 if (fPadIJ[0][j] == fPadIJ[0][i]) continue; // same cathode
799 if (TMath::Abs(fXyq[0][j]-fXyq[0][i]) > 1.e-4) continue;
800 if (TMath::Abs(fXyq[1][j]-fXyq[1][i]) > 1.e-4) continue;
801 fXyq[2][fnPads[0]] += fXyq[2][j];
804 if (cath) fXyq[5][fnPads[0]] = fXyq[5][j]; // save digit number for cath 0
807 // Flag that the digit from the other cathode
808 if (cath && div == 1) fXyq[5][fnPads[0]] = -fXyq[5][i] - 1;
809 // If low pad charge take the other equal to 0
810 if (div == 1 && fXyq[2][fnPads[0]] < fResponse->ZeroSuppression() + 1.5*3) div = 2;
811 fXyq[2][fnPads[0]] /= div;
812 fXyq[0][fnPads[0]] = fXyq[0][i];
813 fXyq[1][fnPads[0]] = fXyq[1][i];
814 fPadIJ[0][fnPads[0]++] = 0;
818 // Check if one-cathode precluster
819 i1 = fnPads[0]!=0 ? 0 : 1;
820 i2 = fnPads[1]!=0 ? 1 : 0;
822 if (i1 != i2) { // two-cathode
824 Int_t *flags = new Int_t[npad];
825 for (Int_t i=0; i<npad; i++) { flags[i] = 0; }
827 // Check pad overlaps
828 for (Int_t i=0; i<npad; i++) {
829 if (fPadIJ[0][i] != i1) continue;
830 xy1[0] = fXyq[0][i] - fXyq[3][i];
831 xy1[1] = fXyq[0][i] + fXyq[3][i];
832 xy1[2] = fXyq[1][i] - fXyq[4][i];
833 xy1[3] = fXyq[1][i] + fXyq[4][i];
834 for (Int_t j=0; j<npad; j++) {
835 if (fPadIJ[0][j] != i2) continue;
836 if (!Overlap(xy1, j, xy12, 0)) continue;
837 flags[i] = flags[j] = 1; // mark overlapped pads
841 // Check if all pads overlap
843 for (Int_t i=0; i<npad; i++) {
844 if (flags[i]) continue;
846 if (fDebug) cout << i << " " << fPadIJ[0][i] << endl;
848 if (fDebug && nFlags) cout << " nFlags = " << nFlags << endl;
849 //if (nFlags > 2 || (Float_t)nFlags / npad > 0.2) { // why 2 ??? - empirical choice
851 for (Int_t i=0; i<npad; i++) {
852 if (flags[i]) continue;
853 digit = TMath::Nint (fXyq[5][i]);
855 fUsed[cath][digit] = kFALSE; // release pad
861 // Check correlations of cathode charges
862 if (fnPads[0] && fnPads[1]) { // two-cathode
864 Int_t over[2] = {1, 1};
865 for (Int_t i=0; i<npad; i++) {
867 if (fXyq[2][i] > 0) sum[cath] += fXyq[2][i];
868 //AZ if (fXyq[2][i] > fResponse->MaxAdc()-1) over[cath] = 0;
869 if (fXyq[2][i] > fResponse->Saturation()-1) over[cath] = 0;
871 if (fDebug) cout << " Total charge: " << sum[0] << " " << sum[1] << endl;
872 if ((over[0] || over[1]) && TMath::Abs(sum[0]-sum[1])/(sum[0]+sum[1])*2 > 1) { // 3 times difference
873 if (fDebug) cout << " Release " << endl;
875 cath = sum[0]>sum[1] ? 0 : 1;
878 Double_t *dist = new Double_t[npad];
879 for (Int_t i=0; i<npad; i++) {
880 if (fPadIJ[0][i] != cath) continue;
881 if (fXyq[2][i] < cmax) continue;
885 // Arrange pads according to their distance to the max,
886 // normalized to the pad size
887 for (Int_t i=0; i<npad; i++) {
889 if (fPadIJ[0][i] != cath) continue;
890 if (i == imax) continue;
891 if (fXyq[2][i] < 0) continue;
892 dist[i] = (fXyq[0][i]-fXyq[0][imax])*(fXyq[0][i]-fXyq[0][imax])/
893 fXyq[3][imax]/fXyq[3][imax]/4;
894 dist[i] += (fXyq[1][i]-fXyq[1][imax])*(fXyq[1][i]-fXyq[1][imax])/
895 fXyq[4][imax]/fXyq[4][imax]/4;
896 dist[i] = TMath::Sqrt (dist[i]);
898 TMath::Sort(npad, dist, flags, kFALSE); // in increasing order
901 for (Int_t i=0; i<npad; i++) {
903 if (fPadIJ[0][indx] != cath) continue;
904 if (fXyq[2][indx] < 0) continue;
905 if (fXyq[2][indx] <= cmax || TMath::Abs(dist[indx]-xmax)<1.e-3) {
907 if (TMath::Abs(dist[indx]-xmax)<1.e-3)
908 cmax = TMath::Max((Double_t)(fXyq[2][indx]),cmax);
909 else cmax = fXyq[2][indx];
911 digit = TMath::Nint (fXyq[5][indx]);
912 fUsed[cath][digit] = kFALSE;
915 // xmax = dist[i]; // Bug?
917 // Check pad overlaps once more
918 for (Int_t i=0; i<npad; i++) flags[i] = 0;
919 for (Int_t i=0; i<npad; i++) {
920 if (fXyq[2][i] < 0) continue;
921 if (fPadIJ[0][i] != i1) continue;
922 xy1[0] = fXyq[0][i] - fXyq[3][i];
923 xy1[1] = fXyq[0][i] + fXyq[3][i];
924 xy1[2] = fXyq[1][i] - fXyq[4][i];
925 xy1[3] = fXyq[1][i] + fXyq[4][i];
926 for (Int_t j=0; j<npad; j++) {
927 if (fXyq[2][j] < 0) continue;
928 if (fPadIJ[0][j] != i2) continue;
929 if (!Overlap(xy1, j, xy12, 0)) continue;
930 flags[i] = flags[j] = 1; // mark overlapped pads
934 for (Int_t i=0; i<npad; i++) {
935 if (fXyq[2][i] < 0 || flags[i]) continue;
938 if (nFlags == fnPads[0] + fnPads[1]) {
940 for (Int_t i=0; i<npad; i++) {
941 if (fXyq[2][i] < 0 || fPadIJ[0][i] != cath) continue;
949 delete [] dist; dist = 0;
950 } // TMath::Abs(sum[0]-sum[1])...
951 } // if (fnPads[0] && fnPads[1])
952 delete [] flags; flags = 0;
955 if (!sameSize) { nShown[0] += fnPads[0]; nShown[1] += fnPads[1]; }
957 // Move released pads to the right
958 Int_t beg = 0, end = npad-1, padij;
961 if (fXyq[2][beg] > 0) { beg++; continue; }
962 for (Int_t j=end; j>beg; j--) {
963 if (fXyq[2][j] < 0) continue;
965 for (Int_t j1=0; j1<2; j1++) {
966 padij = fPadIJ[j1][beg];
967 fPadIJ[j1][beg] = fPadIJ[j1][j];
968 fPadIJ[j1][j] = padij;
970 for (Int_t j1=0; j1<6; j1++) {
972 fXyq[j1][beg] = fXyq[j1][j];
976 } // for (Int_t j=end;
979 npad = fnPads[0] + fnPads[1];
980 if (npad > 500) { cout << " ***** Too large cluster. Give up. " << npad << endl; return kFALSE; }
981 // Back up charge value
982 for (Int_t j=0; j<npad; j++) fXyq[6][j] = fXyq[2][j];
987 //_____________________________________________________________________________
988 void AliMUONClusterFinderAZ::BuildPixArray()
990 // Build pixel array for MLEM method
992 Int_t nPix=0, i1, i2;
993 Float_t xy1[4], xy12[4];
994 AliMUONPixel *pixPtr=0;
996 Int_t npad = fnPads[0] + fnPads[1];
998 // One cathode is empty
999 i1 = fnPads[0]!=0 ? 0 : 1;
1000 i2 = fnPads[1]!=0 ? 1 : 0;
1002 // Build array of pixels on anode plane
1003 if (i1 == i2) { // one-cathode precluster
1004 for (Int_t j=0; j<npad; j++) {
1005 pixPtr = new AliMUONPixel();
1006 for (Int_t i=0; i<2; i++) {
1007 pixPtr->SetCoord(i, fXyq[i][j]); // pixel coordinates
1008 pixPtr->SetSize(i, fXyq[i+3][j]); // pixel size
1010 pixPtr->SetCharge(fXyq[2][j]); // charge
1011 fPixArray->Add((TObject*)pixPtr);
1014 } else { // two-cathode precluster
1015 for (Int_t i=0; i<npad; i++) {
1016 if (fPadIJ[0][i] != i1) continue;
1017 xy1[0] = fXyq[0][i] - fXyq[3][i];
1018 xy1[1] = fXyq[0][i] + fXyq[3][i];
1019 xy1[2] = fXyq[1][i] - fXyq[4][i];
1020 xy1[3] = fXyq[1][i] + fXyq[4][i];
1021 for (Int_t j=0; j<npad; j++) {
1022 if (fPadIJ[0][j] != i2) continue;
1023 if (!Overlap(xy1, j, xy12, 1)) continue;
1024 pixPtr = new AliMUONPixel();
1025 for (Int_t k=0; k<2; k++) {
1026 pixPtr->SetCoord(k, (xy12[2*k]+xy12[2*k+1])/2); // pixel coordinates
1027 pixPtr->SetSize(k, xy12[2*k+1]-pixPtr->Coord(k)); // size
1029 pixPtr->SetCharge(TMath::Min (fXyq[2][i],fXyq[2][j])); //charge
1030 fPixArray->Add((TObject*)pixPtr);
1031 //cout << nPix << " " << pixPtr->Coord(0) << " " << pixPtr->Size(0) << " " << pixPtr->Coord(1) << " " << pixPtr->Size(1) << " " << pixPtr->Charge() << endl;
1033 } // for (Int_t j=0;
1034 } // for (Int_t i=0;
1037 Float_t wxmin=999, wymin=999;
1038 for (Int_t i=0; i<npad; i++) {
1039 //if (fPadIJ[0][i] == i1) wymin = TMath::Min (wymin,fXyq[4][i]);
1040 //if (fPadIJ[0][i] == i2) wxmin = TMath::Min (wxmin,fXyq[3][i]);
1041 wymin = TMath::Min (wymin,fXyq[4][i]);
1042 wxmin = TMath::Min (wxmin,fXyq[3][i]);
1044 if (fDebug) cout << wxmin << " " << wymin << endl;
1046 // Check if small pixel X-size
1047 AdjustPixel(wxmin, 0);
1048 // Check if small pixel Y-size
1049 AdjustPixel(wymin, 1);
1050 // Check if large pixel size
1051 AdjustPixel(wxmin, wymin);
1053 // Remove discarded pixels
1054 for (Int_t i=0; i<nPix; i++) {
1055 pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1057 if (pixPtr->Charge() < 1) { fPixArray->RemoveAt(i); delete pixPtr; }// discarded pixel
1059 fPixArray->Compress();
1060 nPix = fPixArray->GetEntriesFast();
1063 if (fDebug) cout << nPix << endl;
1064 // Too many pixels - sort and remove pixels with the lowest signal
1066 for (Int_t i=npad; i<nPix; i++) {
1067 pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1069 fPixArray->RemoveAt(i);
1073 } // if (nPix > npad)
1075 // Set pixel charges to the same value (for MLEM)
1076 for (Int_t i=0; i<nPix; i++) {
1077 pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1078 //pixPtr->SetCharge(10);
1079 if (fDebug) cout << i+1 << " " << pixPtr->Coord(0) << " " << pixPtr->Coord(1) << " " << pixPtr->Size(0) << " " << pixPtr->Size(1) << endl;
1083 //_____________________________________________________________________________
1084 void AliMUONClusterFinderAZ::AdjustPixel(Float_t width, Int_t ixy)
1086 // Check if some pixels have small size (adjust if necessary)
1088 AliMUONPixel *pixPtr, *pixPtr1 = 0;
1089 Int_t ixy1 = TMath::Even(ixy);
1090 Int_t nPix = fPixArray->GetEntriesFast();
1092 for (Int_t i=0; i<nPix; i++) {
1093 pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1094 if (pixPtr->Charge() < 1) continue; // discarded pixel
1095 if (pixPtr->Size(ixy)-width < -1.e-4) {
1097 if (fDebug) cout << i << " Small X or Y: " << ixy << " " << pixPtr->Size(ixy) << " " << width << " " << pixPtr->Coord(0) << " " << pixPtr->Coord(1) << endl;
1098 for (Int_t j=i+1; j<nPix; j++) {
1099 pixPtr1 = (AliMUONPixel*) fPixArray->UncheckedAt(j);
1100 if (pixPtr1->Charge() < 1) continue; // discarded pixel
1101 if (TMath::Abs(pixPtr1->Size(ixy)-width) < 1.e-4) continue; // right size
1102 if (TMath::Abs(pixPtr1->Coord(ixy1)-pixPtr->Coord(ixy1)) > 1.e-4) continue; // different rows/columns
1103 if (TMath::Abs(pixPtr1->Coord(ixy)-pixPtr->Coord(ixy)) < 2*width) {
1105 //AZ-problem in slats for new segment. pixPtr->SetCoord(ixy, (pixPtr->Coord(ixy)+pixPtr1->Coord(ixy))/2);
1106 Double_t tmp = pixPtr->Coord(ixy) + pixPtr1->Size(ixy) *
1107 TMath::Sign (1., pixPtr1->Coord(ixy) - pixPtr->Coord(ixy));
1108 pixPtr->SetCoord(ixy, tmp);
1109 pixPtr->SetSize(ixy, width);
1110 pixPtr->SetCharge(TMath::Min (pixPtr->Charge(),pixPtr1->Charge()));
1111 pixPtr1->SetCharge(0);
1115 } // for (Int_t j=i+1;
1116 //if (!pixPtr1) { cout << " I am here!" << endl; pixPtr->SetSize(ixy, width); } // ???
1117 //else if (pixPtr1->Charge() > 0.5 || i == nPix-1) {
1118 if (pixPtr1 || i == nPix-1) {
1119 // edge pixel - just increase its size
1120 if (fDebug) cout << " Edge ..." << endl;
1121 for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
1122 //if (fPadIJ[0][j] != ixy1) continue;
1123 if (TMath::Abs(pixPtr->Coord(ixy1)-fXyq[ixy1][j]) > 1.e-4) continue;
1124 if (pixPtr->Coord(ixy) < fXyq[ixy][j])
1125 //pixPtr->Shift(ixy, -pixPtr->Size(ixy));
1126 pixPtr->Shift(ixy, pixPtr->Size(ixy)-width);
1127 //else pixPtr->Shift(ixy, pixPtr->Size(ixy));
1128 else pixPtr->Shift(ixy, -pixPtr->Size(ixy)+width);
1129 pixPtr->SetSize(ixy, width);
1133 } // if (pixPtr->Size(ixy)-width < -1.e-4)
1134 } // for (Int_t i=0; i<nPix;
1138 //_____________________________________________________________________________
1139 void AliMUONClusterFinderAZ::AdjustPixel(Float_t wxmin, Float_t wymin)
1141 // Check if some pixels have large size (adjust if necessary)
1144 Int_t nPix = fPixArray->GetEntriesFast();
1145 AliMUONPixel *pixPtr, *pixPtr1, pix;
1147 // Check if large pixel size
1148 for (Int_t i=0; i<nPix; i++) {
1149 pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1150 if (pixPtr->Charge() < 1) continue; // discarded pixel
1151 if (pixPtr->Size(0)-wxmin > 1.e-4 || pixPtr->Size(1)-wymin > 1.e-4) {
1152 if (fDebug) cout << " Different " << pixPtr->Size(0) << " " << wxmin << " " << pixPtr->Size(1) << " " << wymin << endl;
1154 nx = TMath::Nint (pix.Size(0)/wxmin);
1155 ny = TMath::Nint (pix.Size(1)/wymin);
1156 pix.Shift(0, -pix.Size(0)-wxmin);
1157 pix.Shift(1, -pix.Size(1)-wymin);
1158 pix.SetSize(0, wxmin);
1159 pix.SetSize(1, wymin);
1160 for (Int_t ii=0; ii<nx; ii++) {
1161 pix.Shift(0, wxmin*2);
1162 for (Int_t jj=0; jj<ny; jj++) {
1163 pix.Shift(1, wymin*2);
1164 pixPtr1 = new AliMUONPixel(pix);
1165 fPixArray->Add((TObject*)pixPtr1);
1168 pixPtr->SetCharge(0);
1170 } // for (Int_t i=0; i<nPix;
1174 //_____________________________________________________________________________
1175 Bool_t AliMUONClusterFinderAZ::MainLoop(Int_t iSimple)
1177 // Repeat MLEM algorithm until pixel size becomes sufficiently small
1182 //Int_t nn, xList[10], yList[10];
1183 Int_t nPix = fPixArray->GetEntriesFast();
1184 AliMUONPixel *pixPtr = 0;
1185 Double_t *coef = 0, *probi = 0;
1186 AddVirtualPad(); // add virtual pads if necessary
1187 Int_t npadTot = fnPads[0] + fnPads[1], npadOK = 0;
1188 for (Int_t i=0; i<npadTot; i++) if (fPadIJ[1][i] == 0) npadOK++;
1192 mlem = (TH2D*) gROOT->FindObject("mlem");
1193 if (mlem) mlem->Delete();
1194 // Calculate coefficients
1195 if (fDebug) cout << " nPix, npadTot, npadOK " << nPix << " " << npadTot << " " << npadOK << endl;
1197 // Calculate coefficients and pixel visibilities
1198 coef = new Double_t [npadTot*nPix];
1199 probi = new Double_t [nPix];
1200 for (Int_t ipix=0; ipix<nPix; ipix++) probi[ipix] = 0;
1201 Int_t indx = 0, indx1 = 0, cath = 0;
1203 for (Int_t j=0; j<npadTot; j++) {
1205 if (fPadIJ[1][j] == 0) {
1206 cath = fPadIJ[0][j];
1207 fSegmentation[cath]->GetPadI(fInput->DetElemId(),fXyq[0][j],fXyq[1][j],fZpad,ix,iy);
1208 fSegmentation[cath]->SetPad(fInput->DetElemId(),ix,iy);
1210 fSegmentation[cath]->Neighbours(fInput->DetElemId(),ix,iy,&nn,xList,yList);
1213 for (Int_t i=0; i<nn; i++) {cout << xList[i] << " " << yList[i] << ", ";}
1219 for (Int_t ipix=0; ipix<nPix; ipix++) {
1220 indx1 = indx + ipix;
1221 if (fPadIJ[1][j] < 0) { coef[indx1] = 0; continue; }
1222 pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(ipix);
1223 fSegmentation[cath]->SetHit(fInput->DetElemId(),pixPtr->Coord(0),pixPtr->Coord(1),fZpad);
1224 coef[indx1] = fResponse->IntXY(fInput->DetElemId(),fSegmentation[cath]);
1225 probi[ipix] += coef[indx1];
1226 } // for (Int_t ipix=0;
1227 } // for (Int_t j=0;
1228 for (Int_t ipix=0; ipix<nPix; ipix++) if (probi[ipix] < 0.01) pixPtr->SetCharge(0); // "invisible" pixel
1231 Mlem(coef, probi, 15);
1233 Double_t xylim[4] = {999, 999, 999, 999};
1234 for (Int_t ipix=0; ipix<nPix; ipix++) {
1235 pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(ipix);
1236 for (Int_t i=0; i<4; i++)
1237 xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
1238 //cout << ipix+1; pixPtr->Print();
1240 for (Int_t i=0; i<4; i++) {
1241 xylim[i] -= pixPtr->Size(i/2); if (fDebug) cout << (i%2 ? -1 : 1)*xylim[i] << " "; }
1242 if (fDebug) cout << endl;
1244 // Adjust histogram to approximately the same limits as for the pads
1245 // (for good presentation)
1249 xypads[0] = fHist[0]->GetXaxis()->GetXmin();
1250 xypads[1] = -fHist[0]->GetXaxis()->GetXmax();
1251 xypads[2] = fHist[0]->GetYaxis()->GetXmin();
1252 xypads[3] = -fHist[0]->GetYaxis()->GetXmax();
1253 for (Int_t i=0; i<4; i++) {
1255 if (xylim[i] < xypads[i]) break;
1256 xylim[i] -= 2*pixPtr->Size(i/2);
1261 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
1262 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
1264 mlem = new TH2D("mlem","mlem",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
1265 for (Int_t ipix=0; ipix<nPix; ipix++) {
1266 pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(ipix);
1267 mlem->Fill(pixPtr->Coord(0),pixPtr->Coord(1),pixPtr->Charge());
1269 //gPad->GetCanvas()->cd(3);
1271 ((TCanvas*)gROOT->FindObject("c2"))->cd();
1274 //mlem->SetFillColor(19);
1275 mlem->Draw("lego1Fb");
1280 // Check if the total charge of pixels is too low
1282 for (Int_t i=0; i<nPix; i++) {
1283 pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1284 qTot += pixPtr->Charge();
1286 //AZif (qTot < 1.e-4 || npadOK < 3 && qTot < 50) {
1287 if (qTot < 1.e-4 || npadOK < 3 && qTot < 7) {
1288 delete [] coef; delete [] probi; coef = 0; probi = 0;
1289 fPixArray->Delete();
1290 for (Int_t i=0; i<npadTot; i++) if (fPadIJ[1][i] == 0) fPadIJ[1][i] = -1;
1294 // Plot data - expectation
1296 Double_t x, y, cont;
1297 for (Int_t j=0; j<npadTot; j++) {
1299 for (Int_t i=0; i<nPix; i++) {
1300 // Caculate expectation
1301 pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1302 sum1 += pixPtr->Charge()*coef[j*nPix+i];
1304 //AZsum1 = TMath::Min (sum1,(Double_t)fResponse->MaxAdc());
1305 sum1 = TMath::Min (sum1,(Double_t)fResponse->Saturation());
1308 cath = fPadIJ[0][j];
1309 Int_t ihist = cath*2;
1310 ix = fHist[ihist]->GetXaxis()->FindBin(x);
1311 iy = fHist[ihist]->GetYaxis()->FindBin(y);
1312 cont = fHist[ihist]->GetCellContent(ix,iy);
1313 if (cont == 0 && fHist[ihist+1]) {
1315 ix = fHist[ihist]->GetXaxis()->FindBin(x);
1316 iy = fHist[ihist]->GetYaxis()->FindBin(y);
1318 fHist[ihist]->SetBinContent(ix,iy,fXyq[2][j]-sum1);
1320 ((TCanvas*)gROOT->FindObject("c1"))->cd(1);
1321 //gPad->SetTheta(55);
1323 //mlem->Draw("lego1");
1325 ((TCanvas*)gROOT->FindObject("c1"))->cd(2);
1330 // Simple cluster - skip further passes thru EM-procedure
1331 fxyMu[0][6] = fxyMu[1][6] = 9999;
1333 delete [] coef; delete [] probi; coef = 0; probi = 0;
1334 fPixArray->Delete();
1338 // Calculate position of the center-of-gravity around the maximum pixel
1340 FindCOG(mlem, xyCOG);
1342 if (TMath::Min(pixPtr->Size(0),pixPtr->Size(1)) < 0.07 && pixPtr->Size(0) > pixPtr->Size(1)) break;
1343 //if (TMath::Min(pixPtr->Size(0),pixPtr->Size(1)) < 0.007 && pixPtr->Size(0) > pixPtr->Size(1)) break;
1344 //if (TMath::Min(pixPtr->Size(0),pixPtr->Size(1)) >= 0.07 || pixPtr->Size(0) < pixPtr->Size(1)) {
1345 // Sort pixels according to the charge
1348 for (Int_t i=0; i<nPix; i++) {
1349 pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1350 cout << i+1; pixPtr->Print();
1353 Double_t pixMin = 0.01*((AliMUONPixel*)fPixArray->UncheckedAt(0))->Charge();
1354 pixMin = TMath::Min (pixMin,50.);
1356 // Decrease pixel size and shift pixels to make them centered at
1358 indx = (pixPtr->Size(0)>pixPtr->Size(1)) ? 0 : 1;
1359 Double_t width = 0, shift[2]={0};
1361 for (Int_t i=0; i<4; i++) xylim[i] = 999;
1362 Int_t nPix1 = nPix; nPix = 0;
1363 for (Int_t ipix=0; ipix<nPix1; ipix++) {
1364 pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(ipix);
1365 if (nPix >= npadOK) { // too many pixels already
1366 fPixArray->RemoveAt(ipix);
1370 if (pixPtr->Charge() < pixMin) { // low charge
1371 fPixArray->RemoveAt(ipix);
1375 for (Int_t i=0; i<2; i++) {
1377 pixPtr->SetCharge(10);
1378 pixPtr->SetSize(indx, pixPtr->Size(indx)/2);
1379 width = -pixPtr->Size(indx);
1380 pixPtr->Shift(indx, width);
1381 // Shift pixel position
1384 for (Int_t j=0; j<2; j++) {
1385 shift[j] = pixPtr->Coord(j) - xyCOG[j];
1386 shift[j] -= ((Int_t)(shift[j]/pixPtr->Size(j)/2))*pixPtr->Size(j)*2;
1388 //cout << ipix << " " << i << " " << shift[0] << " " << shift[1] << endl;
1390 pixPtr->Shift(0, -shift[0]);
1391 pixPtr->Shift(1, -shift[1]);
1393 pixPtr = new AliMUONPixel(*pixPtr);
1394 pixPtr->Shift(indx, -2*width);
1395 fPixArray->Add((TObject*)pixPtr);
1398 for (Int_t i=0; i<4; i++)
1399 xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
1400 } // for (Int_t i=0; i<2;
1402 } // for (Int_t ipix=0;
1404 fPixArray->Compress();
1405 nPix = fPixArray->GetEntriesFast();
1407 // Remove excessive pixels
1408 if (nPix > npadOK) {
1409 for (Int_t ipix=npadOK; ipix<nPix; ipix++) {
1410 pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(ipix);
1411 fPixArray->RemoveAt(ipix);
1415 pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(0);
1416 // add pixels if the maximum is at the limit of pixel area
1417 // start from Y-direction
1419 for (Int_t i=3; i>-1; i--) {
1420 if (nPix < npadOK &&
1421 TMath::Abs((i%2 ? -1 : 1)*xylim[i]-xyCOG[i/2]) < pixPtr->Size(i/2)) {
1422 pixPtr = new AliMUONPixel(*pixPtr);
1423 pixPtr->SetCoord(i/2, xyCOG[i/2]+(i%2 ? 2:-2)*pixPtr->Size(i/2));
1424 j = TMath::Even (i/2);
1425 pixPtr->SetCoord(j, xyCOG[j]);
1426 fPixArray->Add((TObject*)pixPtr);
1432 fPixArray->Compress();
1433 nPix = fPixArray->GetEntriesFast();
1434 delete [] coef; delete [] probi; coef = 0; probi = 0;
1437 // remove pixels with low signal or low visibility
1438 // Cuts are empirical !!!
1439 Double_t thresh = TMath::Max (mlem->GetMaximum()/100.,1.);
1440 thresh = TMath::Min (thresh,50.);
1441 Double_t cmax = -1, charge = 0;
1442 for (Int_t i=0; i<nPix; i++) cmax = TMath::Max (cmax,probi[i]);
1443 //cout << thresh << " " << cmax << " " << cmax*0.9 << endl;
1444 // Mark pixels which should be removed
1445 for (Int_t i=0; i<nPix; i++) {
1446 pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1447 charge = pixPtr->Charge();
1448 if (charge < thresh) pixPtr->SetCharge(-charge);
1449 //else if (cmax > 1.91) {
1450 // if (probi[i] < 1.9) pixPtr->SetCharge(-charge);
1452 //AZ else if (probi[i] < cmax*0.9) pixPtr->SetCharge(-charge);
1453 else if (probi[i] < cmax*0.8) pixPtr->SetCharge(-charge);
1454 //cout << i << " " << pixPtr->Coord(0) << " " << pixPtr->Coord(1) << " " << charge << " " << probi[i] << endl;
1456 // Move charge of removed pixels to their nearest neighbour (to keep total charge the same)
1458 for (Int_t i=0; i<nPix; i++) {
1459 pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1460 charge = pixPtr->Charge();
1461 if (charge > 0) continue;
1462 near = FindNearest(pixPtr);
1463 pixPtr->SetCharge(0);
1464 probi[i] = 0; // make it "invisible"
1465 pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(near);
1466 pixPtr->SetCharge(pixPtr->Charge() + (-charge));
1470 for (Int_t i=0; i<nPix; i++) {
1471 pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1472 ix = mlem->GetXaxis()->FindBin(pixPtr->Coord(0));
1473 iy = mlem->GetYaxis()->FindBin(pixPtr->Coord(1));
1474 mlem->SetBinContent(ix, iy, pixPtr->Charge());
1477 ((TCanvas*)gROOT->FindObject("c2"))->cd();
1480 mlem->Draw("lego1Fb");
1484 fxyMu[0][6] = fxyMu[1][6] = 9999;
1485 // Try to split into clusters
1487 if (mlem->GetSum() < 1) ok = kFALSE;
1488 else Split(mlem, coef);
1489 delete [] coef; delete [] probi; coef = 0; probi = 0;
1490 fPixArray->Delete();
1494 //_____________________________________________________________________________
1495 void AliMUONClusterFinderAZ::Mlem(Double_t *coef, Double_t *probi, Int_t nIter)
1497 // Use MLEM to find pixel charges
1499 Int_t nPix = fPixArray->GetEntriesFast();
1500 Int_t npad = fnPads[0] + fnPads[1];
1501 Double_t *probi1 = new Double_t [nPix];
1502 Double_t probMax = 0;
1504 AliMUONPixel *pixPtr;
1506 for (Int_t ipix=0; ipix<nPix; ipix++) if (probi[ipix] > probMax) probMax = probi[ipix];
1507 for (Int_t iter=0; iter<nIter; iter++) {
1509 for (Int_t ipix=0; ipix<nPix; ipix++) {
1510 // Correct each pixel
1511 if (probi[ipix] < 0.01) continue; // skip "invisible" pixel
1513 //probi1[ipix] = probi[ipix];
1514 probi1[ipix] = probMax;
1515 for (Int_t j=0; j<npad; j++) {
1516 if (fPadIJ[1][j] < 0) continue;
1519 indx = indx1 + ipix;
1520 for (Int_t i=0; i<nPix; i++) {
1521 // Caculate expectation
1522 pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1523 sum1 += pixPtr->Charge()*coef[indx1+i];
1524 } // for (Int_t i=0;
1525 //AZ if (fXyq[2][j] > fResponse->MaxAdc()-1 && sum1 > fResponse->MaxAdc()) { probi1[ipix] -= coef[indx]; continue; } // correct for pad charge overflows
1526 if (fXyq[2][j] > fResponse->Saturation()-1 && sum1 > fResponse->Saturation()) { probi1[ipix] -= coef[indx]; continue; } // correct for pad charge overflows
1527 //cout << sum1 << " " << fXyq[2][j] << " " << coef[j*nPix+ipix] << endl;
1528 if (coef[indx] > 1.e-6) sum += fXyq[2][j]*coef[indx]/sum1;
1529 } // for (Int_t j=0;
1530 pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(ipix);
1531 if (probi1[ipix] > 1.e-6) pixPtr->SetCharge(pixPtr->Charge()*sum/probi1[ipix]);
1532 } // for (Int_t ipix=0;
1533 } // for (Int_t iter=0;
1538 //_____________________________________________________________________________
1539 void AliMUONClusterFinderAZ::FindCOG(TH2D *mlem, Double_t *xyc)
1541 // Calculate position of the center-of-gravity around the maximum pixel
1543 Int_t ixmax, iymax, ix, nsumx=0, nsumy=0, nsum=0;
1544 Int_t i1 = -9, j1 = -9;
1545 mlem->GetMaximumBin(ixmax,iymax,ix);
1546 Int_t nx = mlem->GetNbinsX();
1547 Int_t ny = mlem->GetNbinsY();
1548 Double_t thresh = mlem->GetMaximum()/10;
1549 Double_t x, y, cont, xq=0, yq=0, qq=0;
1551 for (Int_t i=TMath::Max(1,iymax-1); i<=TMath::Min(ny,iymax+1); i++) {
1552 //for (Int_t i=TMath::Max(1,iymax-9); i<=TMath::Min(ny,iymax+9); i++) {
1553 y = mlem->GetYaxis()->GetBinCenter(i);
1554 for (Int_t j=TMath::Max(1,ixmax-1); j<=TMath::Min(nx,ixmax+1); j++) {
1555 //for (Int_t j=TMath::Max(1,ixmax-9); j<=TMath::Min(nx,ixmax+9); j++) {
1556 cont = mlem->GetCellContent(j,i);
1557 if (cont < thresh) continue;
1558 if (i != i1) {i1 = i; nsumy++;}
1559 if (j != j1) {j1 = j; nsumx++;}
1560 x = mlem->GetXaxis()->GetBinCenter(j);
1569 Int_t i2 = 0, j2 = 0;
1572 // one bin in Y - add one more (with the largest signal)
1573 for (Int_t i=TMath::Max(1,iymax-1); i<=TMath::Min(ny,iymax+1); i++) {
1574 if (i == iymax) continue;
1575 for (Int_t j=TMath::Max(1,ixmax-1); j<=TMath::Min(nx,ixmax+1); j++) {
1576 cont = mlem->GetCellContent(j,i);
1579 x = mlem->GetXaxis()->GetBinCenter(j);
1580 y = mlem->GetYaxis()->GetBinCenter(i);
1589 if (i2 != i1) nsumy++;
1590 if (j2 != j1) nsumx++;
1592 } // if (nsumy == 1)
1595 // one bin in X - add one more (with the largest signal)
1597 for (Int_t j=TMath::Max(1,ixmax-1); j<=TMath::Min(nx,ixmax+1); j++) {
1598 if (j == ixmax) continue;
1599 for (Int_t i=TMath::Max(1,iymax-1); i<=TMath::Min(ny,iymax+1); i++) {
1600 cont = mlem->GetCellContent(j,i);
1603 x = mlem->GetXaxis()->GetBinCenter(j);
1604 y = mlem->GetYaxis()->GetBinCenter(i);
1613 if (i2 != i1) nsumy++;
1614 if (j2 != j1) nsumx++;
1616 } // if (nsumx == 1)
1618 xyc[0] = xq/qq; xyc[1] = yq/qq;
1619 if (fDebug) cout << xyc[0] << " " << xyc[1] << " " << qq << " " << nsum << " " << nsumx << " " << nsumy << endl;
1623 //_____________________________________________________________________________
1624 Int_t AliMUONClusterFinderAZ::FindNearest(AliMUONPixel *pixPtr0)
1626 // Find the pixel nearest to the given one
1627 // (algorithm may be not very efficient)
1629 Int_t nPix = fPixArray->GetEntriesFast(), imin = 0;
1630 Double_t rmin = 99999, dx = 0, dy = 0, r = 0;
1631 Double_t xc = pixPtr0->Coord(0), yc = pixPtr0->Coord(1);
1632 AliMUONPixel *pixPtr;
1634 for (Int_t i=0; i<nPix; i++) {
1635 pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1636 if (pixPtr->Charge() < 0.5) continue;
1637 dx = (xc - pixPtr->Coord(0)) / pixPtr->Size(0);
1638 dy = (yc - pixPtr->Coord(1)) / pixPtr->Size(1);
1639 r = dx *dx + dy * dy;
1640 if (r < rmin) { rmin = r; imin = i; }
1645 //_____________________________________________________________________________
1646 void AliMUONClusterFinderAZ::Split(TH2D *mlem, Double_t *coef)
1648 // The main steering function to work with clusters of pixels in anode
1649 // plane (find clusters, decouple them from each other, merge them (if
1650 // necessary), pick up coupled pads, call the fitting function)
1652 Int_t nx = mlem->GetNbinsX();
1653 Int_t ny = mlem->GetNbinsY();
1654 Int_t nPix = fPixArray->GetEntriesFast();
1656 Bool_t *used = new Bool_t[ny*nx];
1658 Int_t nclust = 0, indx, indx1;
1660 for (Int_t i=0; i<ny*nx; i++) used[i] = kFALSE;
1662 TObjArray *clusters[200]={0};
1665 // Find clusters of histogram bins (easier to work in 2-D space)
1666 for (Int_t i=1; i<=ny; i++) {
1667 for (Int_t j=1; j<=nx; j++) {
1668 indx = (i-1)*nx + j - 1;
1669 if (used[indx]) continue;
1670 cont = mlem->GetCellContent(j,i);
1671 if (cont < 0.5) continue;
1672 pix = new TObjArray(20);
1674 pix->Add(BinToPix(mlem,j,i));
1675 AddBin(mlem, i, j, 0, used, pix); // recursive call
1676 if (nclust >= 200) AliFatal(" Too many clusters !!!");
1677 clusters[nclust++] = pix;
1678 } // for (Int_t j=1; j<=nx; j++) {
1679 } // for (Int_t i=1; i<=ny;
1680 if (fDebug) cout << nclust << endl;
1681 delete [] used; used = 0;
1683 // Compute couplings between clusters and clusters to pads
1684 Int_t npad = fnPads[0] + fnPads[1];
1686 // Write out some information for algorithm development
1687 Int_t cath=0, npadx[2]={0}, npady[2]={0};
1688 Double_t xlow[2]={9999,9999}, xhig[2]={-9999,-9999};
1689 Double_t ylow[2]={9999,9999}, yhig[2]={-9999,-9999};
1690 for (Int_t j=0; j<npad; j++) {
1691 if (fXyq[3][j] < 0) continue; // exclude virtual pads
1692 cath = fPadIJ[0][j];
1693 if (fXyq[0][j] < xlow[cath]-0.001) {
1694 if (fXyq[0][j]+fXyq[3][j] <= xlow[cath] && npadx[cath]) npadx[cath]++;
1695 xlow[cath] = fXyq[0][j];
1697 if (fXyq[0][j] > xhig[cath]+0.001) {
1698 if (fXyq[0][j]-fXyq[3][j] >= xhig[cath]) npadx[cath]++;
1699 xhig[cath] = fXyq[0][j];
1701 if (fXyq[1][j] < ylow[cath]-0.001) {
1702 if (fXyq[1][j]+fXyq[4][j] <= ylow[cath] && npady[cath]) npady[cath]++;
1703 ylow[cath] = fXyq[1][j];
1705 if (fXyq[1][j] > yhig[cath]+0.001) {
1706 if (fXyq[1][j]-fXyq[4][j] >= yhig[cath]) npady[cath]++;
1707 yhig[cath] = fXyq[1][j];
1710 //if (lun1) fprintf(lun1," %4d %2d %3d %3d %3d %3d \n",gAlice->GetHeader()->GetEvent(),AliMUONClusterInput::Instance()->Chamber(), npadx[0], npadx[1], npady[0], npady[1]);
1712 // Exclude pads with overflows
1713 for (Int_t j=0; j<npad; j++) {
1714 //AZ if (fXyq[2][j] > fResponse->MaxAdc()-1) fPadIJ[1][j] = -5;
1715 if (fXyq[2][j] > fResponse->Saturation()-1) fPadIJ[1][j] = -5;
1716 else fPadIJ[1][j] = 0;
1719 // Compute couplings of clusters to pads
1720 TMatrixD *aijclupad = new TMatrixD(nclust,npad);
1723 for (Int_t iclust=0; iclust<nclust; iclust++) {
1724 pix = clusters[iclust];
1725 npxclu = pix->GetEntriesFast();
1726 for (Int_t i=0; i<npxclu; i++) {
1727 indx = fPixArray->IndexOf(pix->UncheckedAt(i));
1728 for (Int_t j=0; j<npad; j++) {
1729 if (fPadIJ[1][j] < 0 && fPadIJ[1][j] != -5) continue;
1730 if (coef[j*nPix+indx] < fgkCouplMin) continue;
1731 (*aijclupad)(iclust,j) += coef[j*nPix+indx];
1735 // Compute couplings between clusters
1736 TMatrixD *aijcluclu = new TMatrixD(nclust,nclust);
1738 for (Int_t iclust=0; iclust<nclust; iclust++) {
1739 for (Int_t j=0; j<npad; j++) {
1740 // Exclude overflows
1741 if (fPadIJ[1][j] < 0) continue;
1742 if ((*aijclupad)(iclust,j) < fgkCouplMin) continue;
1743 for (Int_t iclust1=iclust+1; iclust1<nclust; iclust1++) {
1744 if ((*aijclupad)(iclust1,j) < fgkCouplMin) continue;
1745 (*aijcluclu)(iclust,iclust1) +=
1746 TMath::Sqrt ((*aijclupad)(iclust,j)*(*aijclupad)(iclust1,j));
1750 for (Int_t iclust=0; iclust<nclust; iclust++) {
1751 for (Int_t iclust1=iclust+1; iclust1<nclust; iclust1++) {
1752 (*aijcluclu)(iclust1,iclust) = (*aijcluclu)(iclust,iclust1);
1756 if (fDebug && nclust > 1) aijcluclu->Print();
1758 // Find groups of coupled clusters
1759 used = new Bool_t[nclust];
1760 for (Int_t i=0; i<nclust; i++) used[i] = kFALSE;
1761 Int_t *clustNumb = new Int_t[nclust];
1762 Int_t nCoupled, nForFit, minGroup[3], clustFit[3], nfit = 0;
1765 for (Int_t igroup=0; igroup<nclust; igroup++) {
1766 if (used[igroup]) continue;
1767 used[igroup] = kTRUE;
1768 clustNumb[0] = igroup;
1770 // Find group of coupled clusters
1771 AddCluster(igroup, nclust, aijcluclu, used, clustNumb, nCoupled); // recursive
1773 cout << " nCoupled: " << nCoupled << endl;
1774 for (Int_t i=0; i<nCoupled; i++) cout << clustNumb[i] << " "; cout << endl;
1776 fnCoupled = nCoupled;
1778 while (nCoupled > 0) {
1782 for (Int_t i=0; i<nCoupled; i++) clustFit[i] = clustNumb[i];
1784 // Too many coupled clusters to fit - try to decouple them
1785 // Find the lowest coupling of 1, 2, min(3,nLinks/2) pixels with
1786 // all the others in the group
1787 for (Int_t j=0; j<3; j++) minGroup[j] = -1;
1788 Double_t coupl = MinGroupCoupl(nCoupled, clustNumb, aijcluclu, minGroup);
1790 // Flag clusters for fit
1792 while (minGroup[nForFit] >= 0 && nForFit < 3) {
1793 if (fDebug) cout << clustNumb[minGroup[nForFit]] << " ";
1794 clustFit[nForFit] = clustNumb[minGroup[nForFit]];
1795 clustNumb[minGroup[nForFit]] -= 999;
1798 if (fDebug) cout << nForFit << " " << coupl << endl;
1801 // Select pads for fit.
1802 if (SelectPad(nCoupled, nForFit, clustNumb, clustFit, aijclupad) < 3 && nCoupled > 1) {
1804 for (Int_t j=0; j<npad; j++) {
1805 if (TMath::Abs(fPadIJ[1][j]) == 1) fPadIJ[1][j] = 0;
1806 if (TMath::Abs(fPadIJ[1][j]) == -9) fPadIJ[1][j] = -5;
1808 // Merge the failed cluster candidates (with too few pads to fit) with
1809 // the one with the strongest coupling
1810 Merge(nForFit, nCoupled, clustNumb, clustFit, clusters, aijcluclu, aijclupad);
1813 nfit = Fit(nForFit, clustFit, clusters, parOk);
1816 // Subtract the fitted charges from pads with strong coupling and/or
1817 // return pads for further use
1818 UpdatePads(nfit, parOk);
1821 for (Int_t j=0; j<npad; j++) {
1822 if (fPadIJ[1][j] == 1) fPadIJ[1][j] = -1;
1823 if (fPadIJ[1][j] == -9) fPadIJ[1][j] = -5;
1826 // Sort the clusters (move to the right the used ones)
1827 Int_t beg = 0, end = nCoupled - 1;
1829 if (clustNumb[beg] >= 0) { beg++; continue; }
1830 for (Int_t j=end; j>beg; j--) {
1831 if (clustNumb[j] < 0) continue;
1833 indx = clustNumb[beg];
1834 clustNumb[beg] = clustNumb[j];
1835 clustNumb[j] = indx;
1841 nCoupled -= nForFit;
1843 // Remove couplings of used clusters
1844 for (Int_t iclust=nCoupled; iclust<nCoupled+nForFit; iclust++) {
1845 indx = clustNumb[iclust] + 999;
1846 for (Int_t iclust1=0; iclust1<nCoupled; iclust1++) {
1847 indx1 = clustNumb[iclust1];
1848 (*aijcluclu)(indx,indx1) = (*aijcluclu)(indx1,indx) = 0;
1852 // Update the remaining clusters couplings (exclude couplings from
1854 for (Int_t j=0; j<npad; j++) {
1855 if (fPadIJ[1][j] != -1) continue;
1856 for (Int_t iclust=0; iclust<nCoupled; iclust++) {
1857 indx = clustNumb[iclust];
1858 if ((*aijclupad)(indx,j) < fgkCouplMin) continue;
1859 for (Int_t iclust1=iclust+1; iclust1<nCoupled; iclust1++) {
1860 indx1 = clustNumb[iclust1];
1861 if ((*aijclupad)(indx1,j) < fgkCouplMin) continue;
1863 (*aijcluclu)(indx,indx1) -=
1864 TMath::Sqrt ((*aijclupad)(indx,j)*(*aijclupad)(indx1,j));
1865 (*aijcluclu)(indx1,indx) = (*aijcluclu)(indx,indx1);
1869 } // for (Int_t j=0; j<npad;
1870 } // if (nCoupled > 3)
1871 } // while (nCoupled > 0)
1872 } // for (Int_t igroup=0; igroup<nclust;
1874 //delete aij_clu; aij_clu = 0; delete aijclupad; aijclupad = 0;
1875 aijcluclu->Delete(); aijclupad->Delete();
1876 for (Int_t iclust=0; iclust<nclust; iclust++) {
1877 pix = clusters[iclust];
1879 delete pix; pix = 0;
1881 delete [] clustNumb; clustNumb = 0; delete [] used; used = 0;
1884 //_____________________________________________________________________________
1885 void AliMUONClusterFinderAZ::AddBin(TH2D *mlem, Int_t ic, Int_t jc, Int_t mode, Bool_t *used, TObjArray *pix)
1887 // Add a bin to the cluster
1889 Int_t nx = mlem->GetNbinsX();
1890 Int_t ny = mlem->GetNbinsY();
1891 Double_t cont1, cont = mlem->GetCellContent(jc,ic);
1892 AliMUONPixel *pixPtr = 0;
1894 for (Int_t i=TMath::Max(ic-1,1); i<=TMath::Min(ic+1,ny); i++) {
1895 for (Int_t j=TMath::Max(jc-1,1); j<=TMath::Min(jc+1,nx); j++) {
1896 if (i != ic && j != jc) continue;
1897 if (used[(i-1)*nx+j-1]) continue;
1898 cont1 = mlem->GetCellContent(j,i);
1899 if (mode && cont1 > cont) continue;
1900 used[(i-1)*nx+j-1] = kTRUE;
1901 if (cont1 < 0.5) continue;
1902 if (pix) pix->Add(BinToPix(mlem,j,i));
1904 pixPtr = new AliMUONPixel (mlem->GetXaxis()->GetBinCenter(j),
1905 mlem->GetYaxis()->GetBinCenter(i), 0, 0, cont1);
1906 fPixArray->Add((TObject*)pixPtr);
1908 AddBin(mlem, i, j, mode, used, pix); // recursive call
1913 //_____________________________________________________________________________
1914 TObject* AliMUONClusterFinderAZ::BinToPix(TH2D *mlem, Int_t jc, Int_t ic)
1916 // Translate histogram bin to pixel
1918 Double_t yc = mlem->GetYaxis()->GetBinCenter(ic);
1919 Double_t xc = mlem->GetXaxis()->GetBinCenter(jc);
1921 Int_t nPix = fPixArray->GetEntriesFast();
1922 AliMUONPixel *pixPtr;
1924 // Compare pixel and bin positions
1925 for (Int_t i=0; i<nPix; i++) {
1926 pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
1927 if (pixPtr->Charge() < 0.5) continue;
1928 if (TMath::Abs(pixPtr->Coord(0)-xc)<1.e-4 && TMath::Abs(pixPtr->Coord(1)-yc)<1.e-4) return (TObject*) pixPtr;
1930 AliWarning(Form(" Something wrong ??? %f %f ", xc, yc));
1934 //_____________________________________________________________________________
1935 void AliMUONClusterFinderAZ::AddCluster(Int_t ic, Int_t nclust, TMatrixD *aijcluclu, Bool_t *used, Int_t *clustNumb, Int_t &nCoupled)
1937 // Add a cluster to the group of coupled clusters
1939 for (Int_t i=0; i<nclust; i++) {
1940 if (used[i]) continue;
1941 if ((*aijcluclu)(i,ic) < fgkCouplMin) continue;
1943 clustNumb[nCoupled++] = i;
1944 AddCluster(i, nclust, aijcluclu, used, clustNumb, nCoupled);
1948 //_____________________________________________________________________________
1949 Double_t AliMUONClusterFinderAZ::MinGroupCoupl(Int_t nCoupled, Int_t *clustNumb, TMatrixD *aijcluclu, Int_t *minGroup)
1951 // Find group of clusters with minimum coupling to all the others
1953 Int_t i123max = TMath::Min(3,nCoupled/2);
1954 Int_t indx, indx1, indx2, indx3, nTot = 0;
1955 Double_t *coupl1 = 0, *coupl2 = 0, *coupl3 = 0;
1957 for (Int_t i123=1; i123<=i123max; i123++) {
1960 coupl1 = new Double_t [nCoupled];
1961 for (Int_t i=0; i<nCoupled; i++) coupl1[i] = 0;
1963 else if (i123 == 2) {
1964 nTot = nCoupled*nCoupled;
1965 coupl2 = new Double_t [nTot];
1966 for (Int_t i=0; i<nTot; i++) coupl2[i] = 9999;
1968 nTot = nTot*nCoupled;
1969 coupl3 = new Double_t [nTot];
1970 for (Int_t i=0; i<nTot; i++) coupl3[i] = 9999;
1973 for (Int_t i=0; i<nCoupled; i++) {
1974 indx1 = clustNumb[i];
1975 for (Int_t j=i+1; j<nCoupled; j++) {
1976 indx2 = clustNumb[j];
1978 coupl1[i] += (*aijcluclu)(indx1,indx2);
1979 coupl1[j] += (*aijcluclu)(indx1,indx2);
1981 else if (i123 == 2) {
1982 indx = i*nCoupled + j;
1983 coupl2[indx] = coupl1[i] + coupl1[j];
1984 coupl2[indx] -= 2 * ((*aijcluclu)(indx1,indx2));
1986 for (Int_t k=j+1; k<nCoupled; k++) {
1987 indx3 = clustNumb[k];
1988 indx = i*nCoupled*nCoupled + j*nCoupled + k;
1989 coupl3[indx] = coupl2[i*nCoupled+j] + coupl1[k];
1990 coupl3[indx] -= 2 * ((*aijcluclu)(indx1,indx3)+(*aijcluclu)(indx2,indx3));
1993 } // for (Int_t j=i+1;
1994 } // for (Int_t i=0;
1995 } // for (Int_t i123=1;
1997 // Find minimum coupling
1998 Double_t couplMin = 9999;
2001 for (Int_t i123=1; i123<=i123max; i123++) {
2003 locMin = TMath::LocMin(nCoupled, coupl1);
2004 couplMin = coupl1[locMin];
2005 minGroup[0] = locMin;
2006 delete [] coupl1; coupl1 = 0;
2008 else if (i123 == 2) {
2009 locMin = TMath::LocMin(nCoupled*nCoupled, coupl2);
2010 if (coupl2[locMin] < couplMin) {
2011 couplMin = coupl2[locMin];
2012 minGroup[0] = locMin/nCoupled;
2013 minGroup[1] = locMin%nCoupled;
2015 delete [] coupl2; coupl2 = 0;
2017 locMin = TMath::LocMin(nTot, coupl3);
2018 if (coupl3[locMin] < couplMin) {
2019 couplMin = coupl3[locMin];
2020 minGroup[0] = locMin/nCoupled/nCoupled;
2021 minGroup[1] = locMin%(nCoupled*nCoupled)/nCoupled;
2022 minGroup[2] = locMin%nCoupled;
2024 delete [] coupl3; coupl3 = 0;
2026 } // for (Int_t i123=1;
2030 //_____________________________________________________________________________
2031 Int_t AliMUONClusterFinderAZ::SelectPad(Int_t nCoupled, Int_t nForFit, Int_t *clustNumb, Int_t *clustFit, TMatrixD *aijclupad)
2033 // Select pads for fit. If too many coupled clusters, find pads giving
2034 // the strongest coupling with the rest of clusters and exclude them from the fit.
2036 Int_t npad = fnPads[0] + fnPads[1];
2037 Double_t *padpix = 0;
2040 padpix = new Double_t[npad];
2041 for (Int_t i=0; i<npad; i++) padpix[i] = 0;
2044 Int_t nOK = 0, indx, indx1;
2045 for (Int_t iclust=0; iclust<nForFit; iclust++) {
2046 indx = clustFit[iclust];
2047 for (Int_t j=0; j<npad; j++) {
2048 if ((*aijclupad)(indx,j) < fgkCouplMin) continue;
2049 if (fPadIJ[1][j] == -5) fPadIJ[1][j] = -9; // flag overflow
2050 if (fPadIJ[1][j] < 0) continue; // exclude overflows and used pads
2051 if (!fPadIJ[1][j]) { fPadIJ[1][j] = 1; nOK++; } // pad to be used in fit
2053 // Check other clusters
2054 for (Int_t iclust1=0; iclust1<nCoupled; iclust1++) {
2055 indx1 = clustNumb[iclust1];
2056 if (indx1 < 0) continue;
2057 if ((*aijclupad)(indx1,j) < fgkCouplMin) continue;
2058 padpix[j] += (*aijclupad)(indx1,j);
2060 } // if (nCoupled > 3)
2061 } // for (Int_t j=0; j<npad;
2062 } // for (Int_t iclust=0; iclust<nForFit
2063 if (nCoupled < 4) return nOK;
2066 for (Int_t j=0; j<npad; j++) {
2067 if (padpix[j] < fgkCouplMin) continue;
2068 if (fDebug) cout << j << " " << padpix[j] << " " << fXyq[0][j] << " " << fXyq[1][j] << endl;
2070 fPadIJ[1][j] = -1; // exclude pads with strong coupling to the other clusters
2073 delete [] padpix; padpix = 0;
2077 //_____________________________________________________________________________
2078 void AliMUONClusterFinderAZ::Merge(Int_t nForFit, Int_t nCoupled, Int_t *clustNumb, Int_t *clustFit, TObjArray **clusters, TMatrixD *aijcluclu, TMatrixD *aijclupad)
2080 // Merge the group of clusters with the one having the strongest coupling with them
2082 Int_t indx, indx1, npxclu, npxclu1, imax=0;
2083 TObjArray *pix, *pix1;
2086 for (Int_t icl=0; icl<nForFit; icl++) {
2087 indx = clustFit[icl];
2088 pix = clusters[indx];
2089 npxclu = pix->GetEntriesFast();
2091 for (Int_t icl1=0; icl1<nCoupled; icl1++) {
2092 indx1 = clustNumb[icl1];
2093 if (indx1 < 0) continue;
2094 if ((*aijcluclu)(indx,indx1) > couplMax) {
2095 couplMax = (*aijcluclu)(indx,indx1);
2098 } // for (Int_t icl1=0;
2099 /*if (couplMax < fgkCouplMin) {
2100 cout << " Oops " << couplMax << endl;
2102 cout << icl << " " << indx << " " << npxclu << " " << nLinks << endl;
2106 pix1 = clusters[imax];
2107 npxclu1 = pix1->GetEntriesFast();
2109 for (Int_t i=0; i<npxclu; i++) { pix1->Add(pix->UncheckedAt(i)); pix->RemoveAt(i); }
2110 if (fDebug) cout << " New number of pixels: " << npxclu1 << " " << pix1->GetEntriesFast() << endl;
2111 //Add cluster-to-cluster couplings
2112 //aijcluclu->Print();
2113 for (Int_t icl1=0; icl1<nCoupled; icl1++) {
2114 indx1 = clustNumb[icl1];
2115 if (indx1 < 0 || indx1 == imax) continue;
2116 (*aijcluclu)(indx1,imax) += (*aijcluclu)(indx,indx1);
2117 (*aijcluclu)(imax,indx1) = (*aijcluclu)(indx1,imax);
2119 (*aijcluclu)(indx,imax) = (*aijcluclu)(imax,indx) = 0;
2120 //aijcluclu->Print();
2121 //Add cluster-to-pad couplings
2122 for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
2123 if (fPadIJ[1][j] < 0 && fPadIJ[1][j] != -5) continue; // exclude used pads
2124 (*aijclupad)(imax,j) += (*aijclupad)(indx,j);
2125 (*aijclupad)(indx,j) = 0;
2127 } // for (Int_t icl=0; icl<nForFit;
2130 //_____________________________________________________________________________
2131 Int_t AliMUONClusterFinderAZ::Fit(Int_t nfit, Int_t *clustFit, TObjArray **clusters, Double_t *parOk)
2133 // Find selected clusters to selected pad charges
2135 TH2D *mlem = (TH2D*) gROOT->FindObject("mlem");
2136 //Int_t nx = mlem->GetNbinsX();
2137 //Int_t ny = mlem->GetNbinsY();
2138 Double_t xmin = mlem->GetXaxis()->GetXmin() - mlem->GetXaxis()->GetBinWidth(1);
2139 Double_t xmax = mlem->GetXaxis()->GetXmax() + mlem->GetXaxis()->GetBinWidth(1);
2140 Double_t ymin = mlem->GetYaxis()->GetXmin() - mlem->GetYaxis()->GetBinWidth(1);
2141 Double_t ymax = mlem->GetYaxis()->GetXmax() + mlem->GetYaxis()->GetBinWidth(1);
2142 //Double_t qmin = 0, qmax = 1;
2143 Double_t step[3]={0.01,0.002,0.02}, xPad = 0, yPad = 99999;
2144 Double_t qPad[2] = {0}, xyqPad[2] = {0};
2146 // Number of pads to use and number of virtual pads
2147 Int_t npads = 0, nVirtual = 0;
2148 for (Int_t i=0; i<fnPads[0]+fnPads[1]; i++) {
2149 if (fPadIJ[1][i] == -9 || fPadIJ[1][i] == 1) {
2150 if (fPadIJ[0][i]) xyqPad[1] += fXyq[0][i] * fXyq[2][i];
2151 else xyqPad[0] += fXyq[1][i] * fXyq[2][i];
2152 qPad[fPadIJ[0][i]] += fXyq[2][i];
2154 if (fXyq[3][i] < 0) nVirtual++;
2155 if (fPadIJ[1][i] != 1) continue;
2156 if (fXyq[3][i] > 0) npads++;
2157 if (yPad > 9999) { xPad = fXyq[0][i]; yPad = fXyq[1][i]; }
2158 //if (fPadIJ[0][i]) xPad = fXyq[0][i];
2161 for (Int_t i=0; i<nfit; i++) {cout << i+1 << " " << clustFit[i] << " ";}
2162 cout << nfit << endl;
2163 cout << " Number of pads to fit: " << npads << endl;
2167 if (npads < 2) return 0;
2169 Int_t digit = 0, nfit0 = nfit;
2170 AliMUONDigit *mdig = 0;
2171 Int_t tracks[3] = {-1, -1, -1};
2172 for (Int_t cath=0; cath<2; cath++) {
2173 for (Int_t i=0; i<fnPads[0]+fnPads[1]; i++) {
2174 if (fPadIJ[0][i] != cath) continue;
2175 if (fPadIJ[1][i] != 1) continue;
2176 if (fXyq[3][i] < 0) continue; // exclude virtual pads
2177 digit = TMath::Nint (fXyq[5][i]);
2178 if (digit >= 0) mdig = fInput->Digit(cath,digit);
2179 else mdig = fInput->Digit(TMath::Even(cath),-digit-1);
2180 //if (!mdig) mdig = fInput->Digit(TMath::Even(cath),digit);
2181 if (!mdig) continue; // protection for cluster display
2182 if (mdig->Hit() >= 0) {
2183 if (tracks[0] < 0) {
2184 tracks[0] = mdig->Hit();
2185 tracks[1] = mdig->Track(0);
2186 } else if (mdig->Track(0) < tracks[1]) {
2187 tracks[0] = mdig->Hit();
2188 tracks[1] = mdig->Track(0);
2191 //AZif (mdig->Track(1)) {
2192 if (mdig->Track(1) >= 0 && mdig->Track(1) != tracks[1]) {
2193 if (tracks[2] < 0) tracks[2] = mdig->Track(1);
2194 else tracks[2] = TMath::Min (tracks[2], mdig->Track(1));
2197 //cout << mdig->Hit() << " " << mdig->Track(0) << " " << mdig->Track(1) <<endl;
2198 } // for (Int_t i=0;
2199 } // for (Int_t cath=0;
2200 //cout << tracks[0] << " " << tracks[1] << " " << tracks[2] <<endl;
2202 // Get number of pads in X and Y
2203 Int_t nInX = 0, nInY;
2204 PadsInXandY(nInX, nInY);
2206 // Take cluster maxima as fitting seeds
2208 AliMUONPixel *pixPtr;
2210 Double_t cont, cmax = 0, xseed = 0, yseed = 0, errOk[8], qq = 0;
2211 Double_t xyseed[3][2], qseed[3], xy_Cand[3][2] = {{0},{0}}, sig_Cand[3][2] = {{0},{0}};
2213 for (Int_t ifit=1; ifit<=nfit; ifit++) {
2215 pix = clusters[clustFit[ifit-1]];
2216 npxclu = pix->GetEntriesFast();
2218 for (Int_t clu=0; clu<npxclu; clu++) {
2219 pixPtr = (AliMUONPixel*) pix->UncheckedAt(clu);
2220 cont = pixPtr->Charge();
2224 xseed = pixPtr->Coord(0);
2225 yseed = pixPtr->Coord(1);
2229 xy_Cand[ifit-1][0] += pixPtr->Coord(0) * cont;
2230 xy_Cand[ifit-1][1] += pixPtr->Coord(1) * cont;
2231 sig_Cand[ifit-1][0] += pixPtr->Coord(0) * pixPtr->Coord(0) * cont;
2232 sig_Cand[ifit-1][1] += pixPtr->Coord(1) * pixPtr->Coord(1) * cont;
2234 xy_Cand[0][0] += pixPtr->Coord(0) * cont;
2235 xy_Cand[0][1] += pixPtr->Coord(1) * cont;
2236 sig_Cand[0][0] += pixPtr->Coord(0) * pixPtr->Coord(0) * cont;
2237 sig_Cand[0][1] += pixPtr->Coord(1) * pixPtr->Coord(1) * cont;
2239 xyseed[ifit-1][0] = xseed;
2240 xyseed[ifit-1][1] = yseed;
2241 qseed[ifit-1] = cmax;
2243 xy_Cand[ifit-1][0] /= qq; // <x>
2244 xy_Cand[ifit-1][1] /= qq; // <y>
2245 sig_Cand[ifit-1][0] = sig_Cand[ifit-1][0]/qq - xy_Cand[ifit-1][0]*xy_Cand[ifit-1][0]; // <x^2> - <x>^2
2246 sig_Cand[ifit-1][0] = sig_Cand[ifit-1][0] > 0 ? TMath::Sqrt (sig_Cand[ifit-1][0]) : 0;
2247 sig_Cand[ifit-1][1] = sig_Cand[ifit-1][1]/qq - xy_Cand[ifit-1][1]*xy_Cand[ifit-1][1]; // <y^2> - <y>^2
2248 sig_Cand[ifit-1][1] = sig_Cand[ifit-1][1] > 0 ? TMath::Sqrt (sig_Cand[ifit-1][1]) : 0;
2249 cout << xy_Cand[ifit-1][0] << " " << xy_Cand[ifit-1][1] << " " << sig_Cand[ifit-1][0] << " " << sig_Cand[ifit-1][1] << endl;
2251 } // for (Int_t ifit=1;
2253 xy_Cand[0][0] /= qq; // <x>
2254 xy_Cand[0][1] /= qq; // <y>
2255 sig_Cand[0][0] = sig_Cand[0][0]/qq - xy_Cand[0][0]*xy_Cand[0][0]; // <x^2> - <x>^2
2256 sig_Cand[0][0] = sig_Cand[0][0] > 0 ? TMath::Sqrt (sig_Cand[0][0]) : 0;
2257 sig_Cand[0][1] = sig_Cand[0][1]/qq - xy_Cand[0][1]*xy_Cand[0][1]; // <y^2> - <y>^2
2258 sig_Cand[0][1] = sig_Cand[0][1] > 0 ? TMath::Sqrt (sig_Cand[0][1]) : 0;
2259 if (fDebug) cout << xy_Cand[0][0] << " " << xy_Cand[0][1] << " " << sig_Cand[0][0] << " " << sig_Cand[0][1] << endl;
2261 Int_t nDof, maxSeed[3];
2262 Double_t fmin, chi2o = 9999, chi2n;
2264 // Try to fit with one-track hypothesis, then 2-track. If chi2/dof is
2265 // lower, try 3-track (if number of pads is sufficient).
2267 TMath::Sort(nfit, qseed, maxSeed, kTRUE); // in decreasing order
2268 nfit = TMath::Min (nfit, (npads + 1) / 3);
2270 if (nInX < 3 && nInY < 3 || nInX == 3 && nInY < 3 || nInX < 3 && nInY == 3) nfit = 1; // not enough pads in each direction
2272 //if (nfit > 1) nfit --;
2273 // One pad per direction
2274 //if (nInX == 1) { step[0] /= 1; xyseed[0][0] = xPad; }
2275 //if (nInY == 1) { step[1] /= 1; xyseed[0][1] = yPad; }
2277 Double_t *gin = 0, func0, func1, param[8], param0[2][8], deriv[2][8], step0[8];
2278 Double_t shift[8], stepMax, derMax, parmin[8], parmax[8], func2[2], shift0;
2279 Double_t delta[8], scMax, dder[8], estim, shiftSave = 0;
2280 Int_t min, max, nCall = 0, memory[8] = {0}, nLoop, idMax = 0, iestMax = 0, nFail;
2281 Double_t rad, dist[3] = {0};
2283 for (Int_t iseed=0; iseed<nfit; iseed++) {
2285 if (iseed) { for (Int_t j=0; j<fNpar; j++) param[j] = parOk[j]; } // for bounded params
2286 for (Int_t j=0; j<3; j++) step0[fNpar+j] = shift[fNpar+j] = step[j];
2287 param[fNpar] = xyseed[maxSeed[iseed]][0];
2288 parmin[fNpar] = xmin;
2289 parmax[fNpar++] = xmax;
2290 param[fNpar] = xyseed[maxSeed[iseed]][1];
2291 parmin[fNpar] = ymin;
2292 parmax[fNpar++] = ymax;
2294 param[fNpar] = fNpar == 4 ? 0.5 : 0.3;
2296 parmax[fNpar++] = 1;
2299 // Try new algorithm
2300 min = nLoop = 1; stepMax = func2[1] = derMax = 999999; nFail = 0;
2304 Fcn1(fNpar, gin, func0, param, 1); nCall++;
2305 //cout << " Func: " << func0 << endl;
2308 for (Int_t j=0; j<fNpar; j++) {
2309 param0[max][j] = param[j];
2310 delta[j] = step0[j];
2311 param[j] += delta[j] / 10;
2312 if (j > 0) param[j-1] -= delta[j-1] / 10;
2313 Fcn1(fNpar, gin, func1, param, 1); nCall++;
2314 deriv[max][j] = (func1 - func0) / delta[j] * 10; // first derivative
2315 //cout << j << " " << deriv[max][j] << endl;
2316 dder[j] = param0[0][j] != param0[1][j] ? (deriv[0][j] - deriv[1][j]) /
2317 (param0[0][j] - param0[1][j]) : 0; // second derivative
2319 param[fNpar-1] -= delta[fNpar-1] / 10;
2320 if (nCall > 2000) break;
2322 min = func2[0] < func2[1] ? 0 : 1;
2323 nFail = min == max ? 0 : nFail + 1;
2325 stepMax = derMax = estim = 0;
2326 for (Int_t j=0; j<fNpar; j++) {
2327 // Estimated distance to minimum
2329 if (nLoop == 1) shift[j] = TMath::Sign (step0[j], -deriv[max][j]); // first step
2330 else if (TMath::Abs(deriv[0][j]) < 1.e-3 && TMath::Abs(deriv[1][j]) < 1.e-3) shift[j] = 0;
2331 else if (deriv[min][j]*deriv[!min][j] > 0 && TMath::Abs(deriv[min][j]) > TMath::Abs(deriv[!min][j])
2332 //|| TMath::Abs(deriv[0][j]-deriv[1][j]) < 1.e-3) {
2333 || TMath::Abs(deriv[0][j]-deriv[1][j]) < 1.e-3 || TMath::Abs(dder[j]) < 1.e-6) {
2334 shift[j] = -TMath::Sign (shift[j], (func2[0]-func2[1]) * (param0[0][j]-param0[1][j]));
2336 if (memory[j] > 1) { shift[j] *= 2; } //cout << " Memory " << memory[j] << " " << shift[j] << endl; }
2340 shift[j] = dder[j] != 0 ? -deriv[min][j] / dder[j] : 0;
2343 if (TMath::Abs(shift[j])/step0[j] > estim) {
2344 estim = TMath::Abs(shift[j])/step0[j];
2349 if (TMath::Abs(shift[j])/step0[j] > 10) shift[j] = TMath::Sign(10.,shift[j]) * step0[j]; //
2351 // Failed to improve minimum
2354 param[j] = param0[min][j];
2355 if (TMath::Abs(shift[j]+shift0) > 0.1*step0[j]) shift[j] = (shift[j] + shift0) / 2;
2356 else shift[j] /= -2;
2360 if (TMath::Abs(shift[j]*deriv[min][j]) > func2[min])
2361 shift[j] = TMath::Sign (func2[min]/deriv[min][j], shift[j]);
2363 // Introduce step relaxation factor
2364 if (memory[j] < 3) {
2365 scMax = 1 + 4 / TMath::Max(nLoop/2.,1.);
2366 if (TMath::Abs(shift0) > 0 && TMath::Abs(shift[j]/shift0) > scMax)
2367 shift[j] = TMath::Sign (shift0*scMax, shift[j]);
2369 param[j] += shift[j];
2370 //AZ Check parameter limits 27-12-2004
2371 if (param[j] < parmin[j]) {
2372 shift[j] = parmin[j] - param[j];
2373 param[j] = parmin[j];
2374 } else if (param[j] > parmax[j]) {
2375 shift[j] = parmax[j] - param[j];
2376 param[j] = parmax[j];
2378 //cout << " xxx " << j << " " << shift[j] << " " << param[j] << endl;
2379 stepMax = TMath::Max (stepMax, TMath::Abs(shift[j]/step0[j]));
2380 if (TMath::Abs(deriv[min][j]) > derMax) {
2382 derMax = TMath::Abs (deriv[min][j]);
2384 } // for (Int_t j=0; j<fNpar;
2385 //cout << max << " " << func2[min] << " " << derMax << " " << stepMax << " " << estim << " " << iestMax << " " << nCall << endl;
2386 if (estim < 1 && derMax < 2 || nLoop > 150) break; // minimum was found
2389 // Check for small step
2390 if (shift[idMax] == 0) { shift[idMax] = step0[idMax]/10; param[idMax] += shift[idMax]; continue; }
2391 if (!memory[idMax] && derMax > 0.5 && nLoop > 10) {
2392 //cout << " ok " << deriv[min][idMax] << " " << deriv[!min][idMax] << " " << dder[idMax]*shift[idMax] << " " << shift[idMax] << endl;
2393 if (dder[idMax] != 0 && TMath::Abs(deriv[min][idMax]/dder[idMax]/shift[idMax]) > 10) {
2394 if (min == max) dder[idMax] = -dder[idMax];
2395 shift[idMax] = -deriv[min][idMax] / dder[idMax] / 10;
2396 param[idMax] += shift[idMax];
2397 stepMax = TMath::Max (stepMax, TMath::Abs(shift[idMax])/step0[idMax]);
2398 //cout << shift[idMax] << " " << param[idMax] << endl;
2399 if (min == max) shiftSave = shift[idMax];
2402 param[idMax] -= shift[idMax];
2403 shift[idMax] = 4 * shiftSave * (gRandom->Rndm(0) - 0.5);
2404 param[idMax] += shift[idMax];
2405 //cout << shift[idMax] << endl;
2411 nDof = npads - fNpar + nVirtual;
2413 chi2n = fmin / nDof;
2414 if (fDebug) cout << " Chi2 " << chi2n << " " << fNpar << endl;
2416 if (chi2n*1.2+1.e-6 > chi2o ) { fNpar -= 3; break; }
2418 // Save parameters and errors
2420 if (nInX == 1 && qPad[1] > 1) {
2421 // One pad per direction
2422 xPad = xyqPad[1] / qPad[1]; // take COG for this case
2423 for (Int_t i=0; i<fNpar; i++) if (i == 0 || i == 2 || i == 5) param0[min][i] = xPad;
2425 if (nInY == 1 && qPad[0] > 1) {
2426 // One pad per direction
2427 yPad = xyqPad[0] / qPad[0]; // take COG for this case
2428 for (Int_t i=0; i<fNpar; i++) if (i == 1 || i == 3 || i == 6) param0[min][i] = yPad;
2433 // Find distance to the nearest neighbour
2434 dist[0] = dist[1] = TMath::Sqrt ((param0[min][0]-param0[min][2])*
2435 (param0[min][0]-param0[min][2])
2436 +(param0[min][1]-param0[min][3])*
2437 (param0[min][1]-param0[min][3]));
2439 dist[2] = TMath::Sqrt ((param0[min][0]-param0[min][5])*
2440 (param0[min][0]-param0[min][5])
2441 +(param0[min][1]-param0[min][6])*
2442 (param0[min][1]-param0[min][6]));
2443 rad = TMath::Sqrt ((param0[min][2]-param0[min][5])*
2444 (param0[min][2]-param0[min][5])
2445 +(param0[min][3]-param0[min][6])*
2446 (param0[min][3]-param0[min][6]));
2447 if (dist[2] < dist[0]) dist[0] = dist[2];
2448 if (rad < dist[1]) dist[1] = rad;
2449 if (rad < dist[2]) dist[2] = rad;
2451 cout << dist[0] << " " << dist[1] << " " << dist[2] << endl;
2452 if (dist[TMath::LocMin(iseed+1,dist)] < 1.) { fNpar -= 3; break; }
2456 for (Int_t i=0; i<fNpar; i++) {
2457 parOk[i] = param0[min][i];
2460 parOk[i] = TMath::Max (parOk[i], parmin[i]);
2461 parOk[i] = TMath::Min (parOk[i], parmax[i]);
2465 if (fmin < 0.1) break; // !!!???
2466 } // for (Int_t iseed=0;
2469 for (Int_t i=0; i<fNpar; i++) {
2470 //if (i == 4 || i == 7) continue;
2471 if (i == 4 || i == 7) {
2472 if (i == 7 || i == 4 && fNpar < 7) cout << parOk[i] << endl;
2473 else cout << parOk[i] * (1-parOk[7]) << endl;
2476 cout << parOk[i] << " " << errOk[i] << endl;
2479 nfit = (fNpar + 1) / 3;
2480 dist[0] = dist[1] = dist[2] = 0;
2483 // Find distance to the nearest neighbour
2484 dist[0] = dist[1] = TMath::Sqrt ((parOk[0]-parOk[2])*
2486 +(parOk[1]-parOk[3])*
2487 (parOk[1]-parOk[3]));
2489 dist[2] = TMath::Sqrt ((parOk[0]-parOk[5])*
2491 +(parOk[1]-parOk[6])*
2492 (parOk[1]-parOk[6]));
2493 rad = TMath::Sqrt ((parOk[2]-parOk[5])*
2495 +(parOk[3]-parOk[6])*
2496 (parOk[3]-parOk[6]));
2497 if (dist[2] < dist[0]) dist[0] = dist[2];
2498 if (rad < dist[1]) dist[1] = rad;
2499 if (rad < dist[2]) dist[2] = rad;
2504 fnPads[1] -= nVirtual;
2507 //for (Int_t j=0; j<nfit; j++) {
2508 for (Int_t j=nfit-1; j>=0; j--) {
2509 indx = j<2 ? j*2 : j*2+1;
2510 if (nfit == 1) coef = 1;
2511 else coef = j==nfit-1 ? parOk[indx+2] : 1-coef;
2512 coef = TMath::Max (coef, 0.);
2513 if (nfit == 3 && j < 2) coef = j==1 ? coef*parOk[indx+2] : coef - parOk[7];
2514 coef = TMath::Max (coef, 0.);
2515 AddRawCluster (parOk[indx], parOk[indx+1], coef*fQtot, errOk[indx], nfit0+10*nfit, tracks,
2516 //sig_Cand[maxSeed[j]][0], sig_Cand[maxSeed[j]][1]);
2517 //sig_Cand[0][0], sig_Cand[0][1], dist[j]);
2518 sig_Cand[0][0], sig_Cand[0][1], dist[TMath::LocMin(nfit,dist)]);
2522 for (Int_t i=0; i<fnMu; i++) {
2524 for (Int_t j=0; j<nfit; j++) {
2525 indx = j<2 ? j*2 : j*2+1;
2526 rad = (fxyMu[i][0]-parOk[indx])*(fxyMu[i][0]-parOk[indx]) +
2527 (fxyMu[i][1]-parOk[indx+1])*(fxyMu[i][1]-parOk[indx+1]);
2532 fxyMu[i][2] = parOk[imax] - fxyMu[i][0];
2533 fxyMu[i][4] = parOk[imax+1] - fxyMu[i][1];
2534 fxyMu[i][3] = errOk[imax];
2535 fxyMu[i][5] = errOk[imax+1];
2542 //_____________________________________________________________________________
2543 void AliMUONClusterFinderAZ::Fcn1(Int_t & /*npar*/, Double_t * /*gin*/, Double_t &f, Double_t *par, Int_t /*iflag*/)
2545 // Fit for one track
2546 //AZ for Muinuit AliMUONClusterFinderAZ& c = *(AliMUONClusterFinderAZ::fgClusterFinder);
2547 AliMUONClusterFinderAZ& c = *this; //AZ
2549 Int_t cath, ix, iy, indx, npads=0;
2550 Double_t charge, delta, coef=0, chi2=0, qTot = 0;
2551 for (Int_t j=0; j<c.fnPads[0]+c.fnPads[1]; j++) {
2552 if (c.fPadIJ[1][j] != 1) continue;
2553 cath = c.fPadIJ[0][j];
2554 if (c.fXyq[3][j] > 0) npads++; // exclude virtual pads
2555 qTot += c.fXyq[2][j];
2556 c.fSegmentation[cath]->GetPadI(fInput->DetElemId(),c.fXyq[0][j],c.fXyq[1][j],c.fZpad,ix,iy);
2557 c.fSegmentation[cath]->SetPad(fInput->DetElemId(),ix,iy);
2559 for (Int_t i=c.fNpar/3; i>=0; i--) { // sum over tracks
2560 indx = i<2 ? 2*i : 2*i+1;
2561 c.fSegmentation[cath]->SetHit(fInput->DetElemId(),par[indx],par[indx+1],c.fZpad);
2562 //charge += c.fResponse->IntXY(fInput->DetElemId(),c.fSegmentation[cath])*par[icl*3+2];
2563 if (c.fNpar == 2) coef = 1;
2564 else coef = i==c.fNpar/3 ? par[indx+2] : 1-coef;
2565 coef = TMath::Max (coef, 0.);
2566 if (c.fNpar == 8 && i < 2) coef = i==1 ? coef*par[indx+2] : coef - par[7];
2567 coef = TMath::Max (coef, 0.);
2568 charge += c.fResponse->IntXY(fInput->DetElemId(),c.fSegmentation[cath])*coef;
2571 //if (c.fXyq[2][j] > c.fResponse->MaxAdc()-1 && charge >
2572 // c.fResponse->MaxAdc()) charge = c.fResponse->MaxAdc();
2573 delta = charge - c.fXyq[2][j];
2575 delta /= c.fXyq[2][j];
2576 //if (cath) delta /= 5; // just for test
2578 } // for (Int_t j=0;
2580 Double_t qAver = qTot/npads; //(c.fnPads[0]+c.fnPads[1]);
2584 //_____________________________________________________________________________
2585 void AliMUONClusterFinderAZ::UpdatePads(Int_t /*nfit*/, Double_t *par)
2587 // Subtract the fitted charges from pads with strong coupling
2589 Int_t cath, ix, iy, indx;
2590 Double_t charge, coef=0;
2591 for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
2592 if (fPadIJ[1][j] != -1) continue;
2594 cath = fPadIJ[0][j];
2595 fSegmentation[cath]->GetPadI(fInput->DetElemId(),fXyq[0][j],fXyq[1][j],fZpad,ix,iy);
2596 fSegmentation[cath]->SetPad(fInput->DetElemId(),ix,iy);
2598 for (Int_t i=fNpar/3; i>=0; i--) { // sum over tracks
2599 indx = i<2 ? 2*i : 2*i+1;
2600 fSegmentation[cath]->SetHit(fInput->DetElemId(),par[indx],par[indx+1],fZpad);
2601 if (fNpar == 2) coef = 1;
2602 else coef = i==fNpar/3 ? par[indx+2] : 1-coef;
2603 coef = TMath::Max (coef, 0.);
2604 if (fNpar == 8 && i < 2) coef = i==1 ? coef*par[indx+2] : coef - par[7];
2605 coef = TMath::Max (coef, 0.);
2606 charge += fResponse->IntXY(fInput->DetElemId(),fSegmentation[cath])*coef;
2609 fXyq[2][j] -= charge;
2610 } // if (fNpar != 0)
2611 if (fXyq[2][j] > fResponse->ZeroSuppression()) fPadIJ[1][j] = 0; // return pad for further using
2612 } // for (Int_t j=0;
2615 //_____________________________________________________________________________
2616 Bool_t AliMUONClusterFinderAZ::TestTrack(Int_t /*t*/) const {
2617 // Test if track was user selected
2620 if (fTrack[0]==-1 || fTrack[1]==-1) {
2622 } else if (t==fTrack[0] || t==fTrack[1]) {
2630 //_____________________________________________________________________________
2631 void AliMUONClusterFinderAZ::AddRawCluster(Double_t x, Double_t y, Double_t qTot, Double_t fmin, Int_t nfit, Int_t *tracks, Double_t /*sigx*/, Double_t /*sigy*/, Double_t /*dist*/)
2634 // Add a raw cluster copy to the list
2636 if (qTot <= 0.501) return;
2637 AliMUONRawCluster cnew;
2638 AliMUON *pMUON=(AliMUON*)gAlice->GetModule("MUON");
2640 //pMUON->AddRawCluster(fInput->Chamber(),c);
2642 Int_t cath, npads[2] = {0}, nover[2] = {0};
2643 for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
2644 cath = fPadIJ[0][j];
2645 // There was an overflow
2646 if (fPadIJ[1][j] == -9) nover[cath]++;
2647 if (fPadIJ[1][j] != 1 && fPadIJ[1][j] != -9) continue;
2648 cnew.SetMultiplicity(cath,cnew.GetMultiplicity(cath)+1);
2649 if (fXyq[2][j] > cnew.GetPeakSignal(cath)) cnew.SetPeakSignal(cath,TMath::Nint (fXyq[2][j]));
2650 //cnew.SetCharge(cath,cnew.GetCharge(cath) + TMath::Nint (fXyq[2][j]));
2651 cnew.SetContrib(npads[cath],cath,fXyq[2][j]);
2652 cnew.SetIndex(npads[cath],cath,TMath::Nint (fXyq[5][j])+10000*fInput->DetElemId());
2656 cnew.SetClusterType(nover[0] + nover[1] * 100);
2657 for (Int_t j=0; j<3; j++) cnew.SetTrack(j,tracks[j]);
2659 for (cath=0; cath<2; cath++) {
2662 cnew.SetZ(cath, fZpad);
2663 cnew.SetCharge(cath, TMath::Nint(qTot));
2664 //cnew.SetPeakSignal(cath,20);
2665 //cnew.SetMultiplicity(cath, 5);
2666 cnew.SetNcluster(cath, nfit);
2667 cnew.SetChi2(cath, fmin); //0.;1
2669 // Evaluate measurement errors
2672 cnew.SetGhost(nfit); //cnew.SetX(1,sigx); cnew.SetY(1,sigy); cnew.SetZ(1,dist);
2673 //cnew.fClusterType=cnew.PhysicsContribution();
2674 //AZ pMUON->GetMUONData()->AddRawCluster(AliMUONClusterInput::Instance()->Chamber(),cnew);
2675 new((*fRawClusters)[fNRawClusters++]) AliMUONRawCluster(cnew); //AZ
2676 if (fDebug) cout << fNRawClusters << " " << AliMUONClusterInput::Instance()->Chamber() << endl;
2680 //_____________________________________________________________________________
2681 Int_t AliMUONClusterFinderAZ::FindLocalMaxima(Int_t *localMax, Double_t *maxVal)
2683 // Find local maxima in pixel space for large preclusters in order to
2684 // try to split them into smaller pieces (to speed up the MLEM procedure)
2686 TH2D *hist = (TH2D*) gROOT->FindObject("anode");
2687 if (hist) hist->Delete();
2689 Double_t xylim[4] = {999, 999, 999, 999};
2690 Int_t nPix = fPixArray->GetEntriesFast();
2691 AliMUONPixel *pixPtr = 0;
2692 for (Int_t ipix=0; ipix<nPix; ipix++) {
2693 pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(ipix);
2694 for (Int_t i=0; i<4; i++)
2695 xylim[i] = TMath::Min (xylim[i], (i%2 ? -1 : 1)*pixPtr->Coord(i/2));
2697 for (Int_t i=0; i<4; i++) xylim[i] -= pixPtr->Size(i/2);
2699 Int_t nx = TMath::Nint ((-xylim[1]-xylim[0])/pixPtr->Size(0)/2);
2700 Int_t ny = TMath::Nint ((-xylim[3]-xylim[2])/pixPtr->Size(1)/2);
2701 hist = new TH2D("anode","anode",nx,xylim[0],-xylim[1],ny,xylim[2],-xylim[3]);
2702 for (Int_t ipix=0; ipix<nPix; ipix++) {
2703 pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(ipix);
2704 hist->Fill(pixPtr->Coord(0), pixPtr->Coord(1), pixPtr->Charge());
2707 ((TCanvas*)gROOT->FindObject("c2"))->cd();
2710 hist->Draw("lego1Fb");
2716 Int_t nMax = 0, indx;
2717 Int_t *isLocalMax = new Int_t[ny*nx];
2718 for (Int_t i=0; i<ny*nx; i++) isLocalMax[i] = 0;
2720 for (Int_t i=1; i<=ny; i++) {
2722 for (Int_t j=1; j<=nx; j++) {
2723 if (hist->GetCellContent(j,i) < 0.5) continue;
2724 //if (isLocalMax[indx+j-1] < 0) continue;
2725 if (isLocalMax[indx+j-1] != 0) continue;
2726 FlagLocalMax(hist, i, j, isLocalMax);
2730 for (Int_t i=1; i<=ny; i++) {
2732 for (Int_t j=1; j<=nx; j++) {
2733 if (isLocalMax[indx+j-1] > 0) {
2734 localMax[nMax] = indx + j - 1;
2735 maxVal[nMax++] = hist->GetCellContent(j,i);
2736 if (nMax > 99) AliFatal(" Too many local maxima !!!");
2740 if (fDebug) cout << " Local max: " << nMax << endl;
2741 delete [] isLocalMax; isLocalMax = 0;
2745 //_____________________________________________________________________________
2746 void AliMUONClusterFinderAZ::FlagLocalMax(TH2D *hist, Int_t i, Int_t j, Int_t *isLocalMax)
2748 // Flag pixels (whether or not local maxima)
2750 Int_t nx = hist->GetNbinsX();
2751 Int_t ny = hist->GetNbinsY();
2752 Int_t cont = TMath::Nint (hist->GetCellContent(j,i));
2755 for (Int_t i1=i-1; i1<i+2; i1++) {
2756 if (i1 < 1 || i1 > ny) continue;
2757 for (Int_t j1=j-1; j1<j+2; j1++) {
2758 if (j1 < 1 || j1 > nx) continue;
2759 if (i == i1 && j == j1) continue;
2760 cont1 = TMath::Nint (hist->GetCellContent(j1,i1));
2761 if (cont < cont1) { isLocalMax[(i-1)*nx+j-1] = -1; return; }
2762 else if (cont > cont1) isLocalMax[(i1-1)*nx+j1-1] = -1;
2763 else { // the same charge
2764 isLocalMax[(i-1)*nx+j-1] = 1;
2765 if (isLocalMax[(i1-1)*nx+j1-1] == 0) {
2766 FlagLocalMax(hist, i1, j1, isLocalMax);
2767 if (isLocalMax[(i1-1)*nx+j1-1] < 0) { isLocalMax[(i-1)*nx+j-1] = -1; return; }
2768 else isLocalMax[(i1-1)*nx+j1-1] = -1;
2773 isLocalMax[(i-1)*nx+j-1] = 1; // local maximum
2776 //_____________________________________________________________________________
2777 void AliMUONClusterFinderAZ::FindCluster(Int_t *localMax, Int_t iMax)
2779 // Find pixel cluster around local maximum #iMax and pick up pads
2780 // overlapping with it
2782 TH2D *hist = (TH2D*) gROOT->FindObject("anode");
2783 Int_t nx = hist->GetNbinsX();
2784 Int_t ny = hist->GetNbinsY();
2785 Int_t ic = localMax[iMax] / nx + 1;
2786 Int_t jc = localMax[iMax] % nx + 1;
2787 Bool_t *used = new Bool_t[ny*nx];
2788 for (Int_t i=0; i<ny*nx; i++) used[i] = kFALSE;
2790 // Drop all pixels from the array - pick up only the ones from the cluster
2791 fPixArray->Delete();
2793 Double_t wx = hist->GetXaxis()->GetBinWidth(1)/2;
2794 Double_t wy = hist->GetYaxis()->GetBinWidth(1)/2;
2795 Double_t yc = hist->GetYaxis()->GetBinCenter(ic);
2796 Double_t xc = hist->GetXaxis()->GetBinCenter(jc);
2797 Double_t cont = hist->GetCellContent(jc,ic);
2798 AliMUONPixel *pixPtr = new AliMUONPixel (xc, yc, wx, wy, cont);
2799 fPixArray->Add((TObject*)pixPtr);
2800 used[(ic-1)*nx+jc-1] = kTRUE;
2801 AddBin(hist, ic, jc, 1, used, (TObjArray*)0); // recursive call
2803 Int_t nPix = fPixArray->GetEntriesFast(), npad = fnPads[0] + fnPads[1];
2804 for (Int_t i=0; i<nPix; i++) {
2805 ((AliMUONPixel*)fPixArray->UncheckedAt(i))->SetSize(0,wx);
2806 ((AliMUONPixel*)fPixArray->UncheckedAt(i))->SetSize(1,wy);
2808 if (fDebug) cout << iMax << " " << nPix << endl;
2810 Float_t xy[4], xy12[4];
2811 // Pick up pads which overlap with found pixels
2812 for (Int_t i=0; i<npad; i++) fPadIJ[1][i] = -1;
2813 for (Int_t i=0; i<nPix; i++) {
2814 pixPtr = (AliMUONPixel*) fPixArray->UncheckedAt(i);
2815 for (Int_t j=0; j<4; j++)
2816 xy[j] = pixPtr->Coord(j/2) + (j%2 ? 1 : -1)*pixPtr->Size(j/2);
2817 for (Int_t j=0; j<npad; j++)
2818 if (Overlap(xy, j, xy12, 0)) fPadIJ[1][j] = 0; // flag for use
2821 delete [] used; used = 0;
2824 //_____________________________________________________________________________
2825 AliMUONClusterFinderAZ&
2826 AliMUONClusterFinderAZ::operator=(const AliMUONClusterFinderAZ& rhs)
2828 // Protected assignement operator
2830 if (this == &rhs) return *this;
2832 AliFatal("Not implemented.");
2837 //_____________________________________________________________________________
2838 void AliMUONClusterFinderAZ::AddVirtualPad()
2840 // Add virtual pad (with small charge) to improve fit for some
2841 // clusters (when pad with max charge is at the extreme of the cluster)
2843 // Get number of pads in X and Y-directions
2844 Int_t nInX = -1, nInY;
2845 PadsInXandY(nInX, nInY);
2849 //nInX = npadx[1] ? npadx[1] : npadx[0];
2850 // Add virtual pads only if number of pads per direction == 2
2851 //if (!npadx[1] && npady[0] != 2 && npadx[0] != 2) return 0; // one-sided
2852 //if (npadx[1] && npady[0] != 2 && npadx[1] != 2) return 0;
2853 if (nInX != 2 && nInY != 2) return;
2855 // Find pads with max charge
2856 Int_t maxpad[2][2] = {{-1, -1}, {-1, -1}}, cath;
2857 Double_t sigmax[2] = {0}, aamax[2] = {0};
2858 for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
2859 if (fPadIJ[1][j] != 0) continue;
2860 cath = fPadIJ[0][j];
2861 if (fXyq[2][j] > sigmax[cath]) {
2862 maxpad[cath][1] = maxpad[cath][0];
2863 aamax[cath] = sigmax[cath];
2864 sigmax[cath] = fXyq[2][j];
2865 maxpad[cath][0] = j;
2868 if (maxpad[0][0] >= 0 && maxpad[0][1] < 0 || maxpad[1][0] >= 0 && maxpad[1][1] < 0) {
2869 for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
2870 if (fPadIJ[1][j] != 0) continue;
2871 cath = fPadIJ[0][j];
2872 if (j == maxpad[cath][0] || j == maxpad[cath][1]) continue;
2873 if (fXyq[2][j] > aamax[cath]) {
2874 aamax[cath] = fXyq[2][j];
2875 maxpad[cath][1] = j;
2879 // Check for mirrors (side X on cathode 0)
2880 Bool_t mirror = kFALSE;
2881 if (maxpad[0][0] >= 0 && maxpad[1][0] >= 0)
2882 mirror = fXyq[3][maxpad[0][0]] < fXyq[4][maxpad[0][0]];
2884 // Find neughbours of pads with max charges
2885 Int_t nn, xList[10], yList[10], ix0, iy0, ix, iy, neighb;
2886 for (cath=0; cath<2; cath++) {
2887 if (!cath && maxpad[0][0] < 0) continue; // one-sided cluster - cathode 1
2888 if (cath && maxpad[1][0] < 0) break; // one-sided cluster - cathode 0
2889 if (maxpad[1][0] >= 0) {
2891 if (!cath && nInY != 2) continue;
2892 //AZ if (cath && nInX != 2) continue;
2893 if (cath && nInX != 2 && (maxpad[0][0] >= 0 || nInY != 2)) continue;
2895 if (!cath && nInX != 2) continue;
2896 if (cath && nInY != 2 && (maxpad[0][0] >= 0 || nInX != 2)) continue;
2900 Int_t iAddX = 0, iAddY = 0, ix1 = 0, iy1 = 0, iMuon = 0, iPad = 0;
2901 if (maxpad[0][0] < 0) iPad = 1;
2903 // This part of code to take care of edge effect (problems in MC)
2904 Float_t spr_x = fResponse->SigmaIntegration()*fResponse->ChargeSpreadX();
2905 Float_t spr_y = fResponse->SigmaIntegration()*fResponse->ChargeSpreadY();
2906 Double_t rmin = 9999, rad2;
2907 Int_t border = 0, iYlow = 0;
2910 for (Int_t i=0; i<2; i++) {
2911 rad2 = (fXyq[0][maxpad[iPad][0]]-fxyMu[i][0]) * (fXyq[0][maxpad[iPad][0]]-fxyMu[i][0]);
2912 rad2 += (fXyq[1][maxpad[iPad][0]]-fxyMu[i][1]) * (fXyq[1][maxpad[iPad][0]]-fxyMu[i][1]);
2913 if (rad2 < rmin) { iMuon = i; rmin = rad2; }
2915 fSegmentation[cath]->FirstPad(fInput->DetElemId(),(Float_t)fxyMu[iMuon][0], (Float_t)fxyMu[iMuon][1], fZpad, spr_x, spr_y);
2916 if (fSegmentation[cath]->Sector(fInput->DetElemId(),fSegmentation[cath]->Ix(),fSegmentation[cath]->Iy()) <= 0) {
2917 fSegmentation[cath]->NextPad(fInput->DetElemId());
2919 iYlow = fSegmentation[cath]->Iy();
2923 for (iPad=0; iPad<2; iPad++) {
2924 if (iPad && !iAddX && !iAddY) break;
2925 if (iPad && fXyq[2][maxpad[cath][1]] / sigmax[cath] < 0.5) break;
2927 Int_t neighbx = 0, neighby = 0;
2928 fSegmentation[cath]->GetPadI(fInput->DetElemId(),fXyq[0][maxpad[cath][iPad]],fXyq[1][maxpad[cath][iPad]],fZpad,ix0,iy0);
2929 fSegmentation[cath]->Neighbours(fInput->DetElemId(),ix0,iy0,&nn,xList,yList);
2930 Float_t zpad; //, xpad, ypad;
2931 for (Int_t j=0; j<nn; j++) {
2933 if (border && yList[j] < iYlow) { xList[j] = yList[j] = 0; continue; }
2934 fSegmentation[cath]->GetPadC(fInput->DetElemId(),xList[j],yList[j],xpad,ypad,zpad);
2935 if (TMath::Abs(xpad) < 1 && TMath::Abs(ypad) < 1)
2936 { xList[j] = yList[j] = 0; continue; } // strange case (something with pad mapping)
2938 if (TMath::Abs(xList[j]-ix0) == 1 || TMath::Abs(xList[j]*ix0) == 1) neighbx++;
2939 if (TMath::Abs(yList[j]-iy0) == 1 || TMath::Abs(yList[j]*iy0) == 1) neighby++;
2942 if (cath) neighb = neighbx;
2943 else neighb = neighby;
2944 if (maxpad[0][0] < 0) neighb += neighby;
2945 else if (maxpad[1][0] < 0) neighb += neighbx;
2947 if (!cath) neighb = neighbx;
2948 else neighb = neighby;
2949 if (maxpad[0][0] < 0) neighb += neighbx;
2950 else if (maxpad[1][0] < 0) neighb += neighby;
2953 for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
2954 if (fPadIJ[0][j] != cath) continue;
2955 fSegmentation[cath]->GetPadI(fInput->DetElemId(),fXyq[0][j],fXyq[1][j],fZpad,ix,iy);
2956 if (iy == iy0 && ix == ix0) continue;
2957 for (Int_t k=0; k<nn; k++) {
2958 if (xList[k] != ix || yList[k] != iy) continue;
2960 if ((!cath || maxpad[0][0] < 0) &&
2961 (TMath::Abs(iy-iy0) == 1 || TMath::Abs(iy*iy0) == 1)) {
2962 xList[k] = yList[k] = 0;
2966 if ((cath || maxpad[1][0] < 0) &&
2967 (TMath::Abs(ix-ix0) == 1 || TMath::Abs(ix*ix0) == 1)) {
2968 xList[k] = yList[k] = 0;
2972 if ((!cath || maxpad[0][0] < 0) &&
2973 (TMath::Abs(ix-ix0) == 1 || TMath::Abs(ix*ix0) == 1)) {
2974 xList[k] = yList[k] = 0;
2978 if ((cath || maxpad[1][0] < 0) &&
2979 (TMath::Abs(iy-iy0) == 1 || TMath::Abs(iy*iy0) == 1)) {
2980 xList[k] = yList[k] = 0;
2985 } // for (Int_t k=0; k<nn;
2987 } // for (Int_t j=0; j<fnPads[0]+fnPads[1];
2988 if (!neighb) continue;
2993 for (Int_t j=0; j<nn; j++) {
2994 if (xList[j] == 0 && yList[j] == 0) continue;
2995 npads = fnPads[0] + fnPads[1];
2996 fPadIJ[0][npads] = cath;
2997 fPadIJ[1][npads] = 0;
3000 if (TMath::Abs(ix-ix0) == 1 || TMath::Abs(ix*ix0) == 1) {
3001 if (iy != iy0) continue; // new segmentation - check
3002 if (nInX != 2) continue; // new
3004 if (!cath && maxpad[1][0] >= 0) continue;
3005 //if (maxpad[1][0] < 0 && nInX != 2) continue;
3007 if (cath && maxpad[0][0] >= 0) continue;
3008 //if (maxpad[0][0] < 0 && nInX != 2) continue;
3010 if (iPad && !iAddX) continue;
3011 fSegmentation[cath]->GetPadC(fInput->DetElemId(),ix,iy,fXyq[0][npads],fXyq[1][npads],zpad);
3012 if (ix1 == ix0) continue;
3013 //if (iPad && ix1 == ix0) continue;
3014 //if (iPad && TMath::Abs(fXyq[0][npads]-fXyq[0][iAddX]) < fXyq[3][iAddX]) continue;
3015 //if (TMath::Abs(fXyq[0][npads]) < 1 && TMath::Abs(fXyq[1][npads]) < 1) continue; // strange case (something with pad mapping)
3016 if (maxpad[1][0] < 0 || mirror && maxpad[0][0] >= 0) {
3017 if (!iPad) fXyq[2][npads] = TMath::Min (sigmax[0]/100, 5.);
3018 else fXyq[2][npads] = TMath::Min (aamax[0]/100, 5.);
3021 if (!iPad) fXyq[2][npads] = TMath::Min (sigmax[1]/100, 5.);
3022 else fXyq[2][npads] = TMath::Min (aamax[1]/100, 5.);
3024 fXyq[2][npads] = TMath::Max (fXyq[2][npads], (float)1);
3025 //fXyq[2][npads] = 1;
3026 //isec = fSegmentation[cath]->Sector(fInput->DetElemId(),ix, iy);
3027 //fXyq[3][npads] = fSegmentation[cath]->Dpx(fInput->DetElemId(),isec)/2;
3028 fXyq[3][npads] = -2; // flag
3031 if (fDebug) cout << " ***** Add virtual pad in X ***** " << fXyq[2][npads]
3032 << " " << fXyq[0][npads] << " " << fXyq[1][npads] << endl;
3036 if (nInY != 2) continue;
3037 if (!mirror && cath && maxpad[0][0] >= 0) continue;
3038 if (mirror && !cath && maxpad[1][0] >= 0) continue;
3039 if (TMath::Abs(iy-iy0) == 1 || TMath::Abs(iy*iy0) == 1) {
3040 if (ix != ix0) continue; // new segmentation - check
3041 if (iPad && !iAddY) continue;
3042 fSegmentation[cath]->GetPadC(fInput->DetElemId(),ix,iy,fXyq[0][npads],fXyq[1][npads],zpad);
3043 if (iy1 == iy0) continue;
3044 //if (iPad && iy1 == iy0) continue;
3045 //if (iPad && TMath::Abs(fXyq[1][npads]-fXyq[1][iAddY]) < fXyq[4][iAddY]) continue;
3046 //if (TMath::Abs(fXyq[0][npads]) < 1 && TMath::Abs(fXyq[1][npads]) < 1) continue; // strange case (something with pad mapping)
3047 if (maxpad[0][0] < 0 || mirror && maxpad[1][0] >= 0) {
3048 if (!iPad) fXyq[2][npads] = TMath::Min (sigmax[1]/20, 5.);
3049 else fXyq[2][npads] = TMath::Min (aamax[1]/20, 5.);
3052 if (!iPad) fXyq[2][npads] = TMath::Min (sigmax[0]/20, 5.);
3053 else fXyq[2][npads] = TMath::Min (aamax[0]/20, 5.);
3055 fXyq[2][npads] = TMath::Max (fXyq[2][npads], (float)1);
3056 //isec = fSegmentation[cath]->Sector(fInput->DetElemId(),ix, iy);
3057 //fXyq[4][npads] = fSegmentation[cath]->Dpy(isec)/2;
3058 fXyq[3][npads] = -2; // flag
3061 if (fDebug) cout << " ***** Add virtual pad in Y ***** " << fXyq[2][npads]
3062 << " " << fXyq[0][npads] << " " << fXyq[1][npads] << endl;
3065 } // for (Int_t j=0; j<nn;
3066 } // for (Int_t iPad=0;
3067 } // for (cath=0; cath<2;
3071 //_____________________________________________________________________________
3072 void AliMUONClusterFinderAZ::PadsInXandY(Int_t &nInX, Int_t &nInY)
3074 // Find number of pads in X and Y-directions (excluding virtual ones and
3077 static Int_t nXsaved = 0, nYsaved = 0;
3078 nXsaved = nYsaved = 0;
3079 //if (nInX >= 0) {nInX = nXsaved; nInY = nYsaved; return; }
3080 Double_t xlow[2] = {9999,9999}, xhig[2] = {-9999,-9999};
3081 Double_t ylow[2] = {9999,9999}, yhig[2] = {-9999,-9999};
3082 Int_t nx, ny, cath, npadx[2] = {0}, npady[2] = {0};
3083 for (Int_t j=0; j<fnPads[0]+fnPads[1]; j++) {
3084 //if (fPadIJ[1][j] != 0) continue;
3085 //if (fXyq[3][j] < 0) continue; // virtual pad
3086 if (nInX < 0 && fPadIJ[1][j] != 0) continue; // before fit
3087 else if (nInX == 0 && fPadIJ[1][j] != 1) continue; // fit - exclude overflows
3088 else if (nInX > 0 && fPadIJ[1][j] != 1 && fPadIJ[1][j] != -9) continue; // exclude non-marked
3089 //AZif (fXyq[2][j] > fResponse->MaxAdc()-1) continue;
3090 if (nInX <= 0 && fXyq[2][j] > fResponse->Saturation()-1) continue; // skip overflows
3091 cath = fPadIJ[0][j];
3093 if (fXyq[0][j] < xlow[cath]-0.001) { xlow[cath] = fXyq[0][j]; nx++; }
3094 if (fXyq[1][j] < ylow[cath]-0.001) { ylow[cath] = fXyq[1][j]; ny++; }
3095 if (fXyq[0][j] > xhig[cath]+0.001) { xhig[cath] = fXyq[0][j]; nx++; }
3096 if (fXyq[1][j] > yhig[cath]+0.001) { yhig[cath] = fXyq[1][j]; ny++; }
3097 if (nx % 2 || !npadx[cath]) npadx[cath]++;
3098 if (ny % 2 || !npady[cath]) npady[cath]++;
3100 //nInY = nYsaved == npady[0] ? npady[0] : npady[1];
3101 //nInX = nXsaved == npadx[1] ? npadx[1] : npadx[0];
3102 nInY = TMath::Max (npady[0], npady[1]);
3103 nInX = TMath::Max (npadx[0], npadx[1]);
3104 //nInY = npady[0] > 0 ? npady[0] : npady[1];
3105 //nInX = npadx[1] > 0 ? npadx[1] : npadx[0];
3108 //_____________________________________________________________________________
3109 void AliMUONClusterFinderAZ::Simple()
3111 // Process simple cluster (small number of pads) without EM-procedure
3113 Int_t nForFit = 1, clustFit[1] = {1}, nfit;
3114 Double_t parOk[3] = {0.};
3115 TObjArray *clusters[1];
3116 clusters[1] = fPixArray;
3117 for (Int_t i=0; i<fnPads[0]+fnPads[1]; i++) fPadIJ[1][i] = 1;
3119 nfit = Fit(nForFit, clustFit, clusters, parOk);
3122 //_____________________________________________________________________________
3123 void AliMUONClusterFinderAZ::Errors(AliMUONRawCluster *clus)
3125 // Correct reconstructed coordinates for some clusters and evaluate errors
3127 Double_t qTot = clus->GetCharge(0), fmin = clus->GetChi2(0);
3128 Double_t xreco = clus->GetX(0), yreco = clus->GetY(0), zreco = clus->GetZ(0);
3129 Double_t sigmax[2] = {0};
3131 Int_t nInX = 1, nInY, maxdig[2] ={-1, -1}, digit, cath1, isec;
3132 PadsInXandY(nInX, nInY);
3134 // Find pad with maximum signal
3135 for (Int_t cath = 0; cath < 2; cath++) {
3136 for (Int_t j = 0; j < clus->GetMultiplicity(cath); j++) {
3138 digit = clus->GetIndex(j, cath);
3139 if (digit < 0) { cath1 = TMath::Even(cath); digit = -digit - 1; } // from the other cathode
3141 if (clus->GetContrib(j,cath) > sigmax[cath1]) {
3142 sigmax[cath1] = clus->GetContrib(j,cath);
3143 maxdig[cath1] = digit;
3148 // Size of pad with maximum signal and reco coordinate distance from the pad center
3149 AliMUONDigit *mdig = 0;
3150 Double_t wx[2], wy[2], dxc[2], dyc[2];
3151 Float_t xpad, ypad, zpad;
3153 for (Int_t cath = 0; cath < 2; cath++) {
3154 if (maxdig[cath] < 0) continue;
3155 mdig = fInput->Digit(cath,maxdig[cath]);
3156 isec = fSegmentation[cath]->Sector(fInput->DetElemId(),mdig->PadX(), mdig->PadY());
3157 wx[cath] = fSegmentation[cath]->Dpx(fInput->DetElemId(),isec);
3158 wy[cath] = fSegmentation[cath]->Dpy(fInput->DetElemId(),isec);
3159 fSegmentation[cath]->GetPadI(fInput->DetElemId(),xreco, yreco, zreco, ix, iy);
3160 isec = fSegmentation[cath]->Sector(fInput->DetElemId(),ix,iy);
3162 fSegmentation[cath]->GetPadC(fInput->DetElemId(), ix, iy, xpad, ypad, zpad);
3163 dxc[cath] = xreco - xpad;
3164 dyc[cath] = yreco - ypad;
3168 // Check if pad with max charge at the edge (number of neughbours)
3169 Int_t nn, xList[10], yList[10], neighbx[2][2] = {{0,0}, {0,0}}, neighby[2][2]= {{0,0}, {0,0}};
3170 for (Int_t cath = 0; cath < 2; cath++) {
3171 if (maxdig[cath] < 0) continue;
3172 mdig = fInput->Digit(cath,maxdig[cath]);
3173 fSegmentation[cath]->Neighbours(fInput->DetElemId(),mdig->PadX(),mdig->PadY(),&nn,xList,yList);
3174 isec = fSegmentation[cath]->Sector(fInput->DetElemId(),mdig->PadX(), mdig->PadY());
3176 Float_t spr_x = fResponse->SigmaIntegration() * fResponse->ChargeSpreadX();
3177 Float_t spr_y = fResponse->SigmaIntegration() * fResponse->ChargeSpreadY();
3178 //fSegmentation[cath]->FirstPad(fInput->DetElemId(),muons[ihit][1], muons[ihit][2], muons[ihit][3], spr_x, spr_y);
3179 fSegmentation[cath]->FirstPad(fInput->DetElemId(),xreco, yreco, zreco, spr_x, spr_y);
3181 if (fSegmentation[cath]->Sector(fInput->DetElemId(),fSegmentation[cath]->Ix(),fSegmentation[cath]->Iy()) <= 0) {
3182 fSegmentation[cath]->NextPad(fInput->DetElemId());
3186 for (Int_t j=0; j<nn; j++) {
3187 if (border && yList[j] < fSegmentation[cath]->Iy()) continue;
3188 fSegmentation[cath]->GetPadC (fInput->DetElemId(), xList[j], yList[j], xpad, ypad, zpad);
3189 //cout << ch << " " << xList[j] << " " << yList[j] << " " << border << " " << x << " " << y << " " << xpad << " " << ypad << endl;
3190 if (TMath::Abs(xpad) < 1 && TMath::Abs(ypad) < 1) continue;
3191 if (xList[j] == mdig->PadX()-1 || mdig->PadX() == 1 &&
3192 xList[j] == -1) neighbx[cath][0] = 1;
3193 else if (xList[j] == mdig->PadX()+1 || mdig->PadX() == -1 &&
3194 xList[j] == 1) neighbx[cath][1] = 1;
3195 if (yList[j] == mdig->PadY()-1 || mdig->PadY() == 1 &&
3196 yList[j] == -1) neighby[cath][0] = 1;
3197 else if (yList[j] == mdig->PadY()+1 || mdig->PadY() == -1 &&
3198 yList[j] == 1) neighby[cath][1] = 1;
3199 } // for (Int_t j=0; j<nn;
3200 if (neighbx[cath][0] && neighbx[cath][1]) neighbx[cath][0] = 0;
3201 else if (neighbx[cath][1]) neighbx[cath][0] = -1;
3202 else neighbx[cath][0] = 1;
3203 if (neighby[cath][0] && neighby[cath][1]) neighby[cath][0] = 0;
3204 else if (neighby[cath][1]) neighby[cath][0] = -1;
3205 else neighby[cath][0] = 1;
3208 Int_t iOver = clus->GetClusterType();
3209 // One-sided cluster
3210 if (!clus->GetMultiplicity(0)) {
3211 neighby[0][0] = neighby[1][0];
3213 if (iOver < 99) iOver += 100 * iOver;
3215 } else if (!clus->GetMultiplicity(1)) {
3216 neighbx[1][0] = neighbx[0][0];
3218 if (iOver < 99) iOver += 100 * iOver;
3222 // Apply corrections and evaluate errors
3223 Double_t errY, errX;
3224 Errors(nInY, nInX, neighby[0][0],neighbx[1][0], fmin, wy[0]*10, wx[1]*10, iOver,
3225 dyc[0], dxc[1], qTot, yreco, xreco, errY, errX);
3226 errY = TMath::Max (errY, 0.01);
3228 //errX = TMath::Max (errX, 0.144);
3229 clus->SetX(0, xreco); clus->SetY(0, yreco);
3230 clus->SetX(1, errX); clus->SetY(1, errY);
3233 //_____________________________________________________________________________
3234 void AliMUONClusterFinderAZ::Errors(Int_t ny, Int_t nx, Int_t iby, Int_t ibx, Double_t fmin,
3235 Double_t wy, Double_t wx, Int_t iover,
3236 Double_t dyc, Double_t /*dxc*/, Double_t qtot,
3237 Double_t &yrec, Double_t &xrec, Double_t &erry, Double_t &errx)
3239 // Correct reconstructed coordinates for some clusters and evaluate errors
3243 Int_t iovery = iover % 100;
3250 yrec += iby * (0.1823+0.2008)/2;
3253 // Find "effective pad width"
3254 Double_t width = 0.218 / (1.31e-4 * TMath::Exp (2.688 * TMath::Log(qtot)) + 1) * 2;
3255 width = TMath::Min (width, 0.4);
3256 erry = width / TMath::Sqrt(12.);
3257 erry = TMath::Max (erry, 0.01293);
3262 /* ---> "Bad" fit */
3265 if (ny == 5) erry = 0.06481;
3272 erry = 0.00417; //0.01010
3275 if (dyc * iby > -0.05) {
3276 Double_t dyc2 = dyc * dyc;
3278 corr = 0.019 - 0.602 * dyc + 8.739 * dyc2 - 44.209 * dyc2 * dyc;
3279 corr = TMath::Min (corr, TMath::Abs(-0.25-dyc));
3284 corr = 0.006 + 0.300 * dyc + 6.147 * dyc2 + 42.039 * dyc2 * dyc;
3285 corr = TMath::Min (corr, 0.25-dyc);
3291 erry = (0.00303 + 0.00296) / 2;
3297 /* ---> Overflows */
3304 } else if (TMath::Abs(wy - 5) < 0.1) {
3305 erry = 0.061; //0.06622
3307 erry = 0.00812; // 0.01073
3313 /* ---> "Good" but very high signal */
3315 if (TMath::Abs(wy - 4) < 0.1) {
3317 } else if (fmin < 0.03 && qtot < 6000) {
3325 /* ---> "Good" clusters */
3327 if (TMath::Abs(wy - 5) < 0.1) {
3328 erry = 0.0011; //0.00304
3329 } else if (qtot < 400.) {
3332 erry = 0.00135; // 0.00358
3334 } else if (ny == 3) {
3335 if (TMath::Abs(wy - 4) < 0.1) {
3336 erry = 35.407 / (1 + TMath::Exp(5.511*TMath::Log(qtot/265.51))) + 11.564;
3337 //erry = 83.512 / (1 + TMath::Exp(3.344*TMath::Log(qtot/211.58))) + 12.260;
3339 erry = 147.03 / (1 + TMath::Exp(1.713*TMath::Log(qtot/73.151))) + 9.575;
3340 //erry = 91.743 / (1 + TMath::Exp(2.332*TMath::Log(qtot/151.67))) + 11.453;
3345 if (TMath::Abs(wy - 4) < 0.1) {
3346 erry = 60.800 / (1 + TMath::Exp(3.305*TMath::Log(qtot/104.53))) + 11.702;
3347 //erry = 73.128 / (1 + TMath::Exp(5.676*TMath::Log(qtot/120.93))) + 17.839;
3349 erry = 117.98 / (1 + TMath::Exp(2.005*TMath::Log(qtot/37.649))) + 21.431;
3350 //erry = 99.066 / (1 + TMath::Exp(4.900*TMath::Log(qtot/107.57))) + 25.315;
3357 /* ---> X-coordinate */
3366 if (TMath::Abs(wx - 6) < 0.1) {
3367 if (qtot < 40) errx = 0.1693;
3368 else errx = 0.06241;
3369 } else if (TMath::Abs(wx - 7.5) < 0.1) {
3370 if (qtot < 40) errx = 0.2173;
3371 else errx = 0.07703;
3372 } else if (TMath::Abs(wx - 10) < 0.1) {
3374 if (qtot < 40) errx = 0.2316;
3377 xrec += (0.2115 + 0.1942) / 2 * ibx;
3383 /* ---> "Bad" fit */
3390 if (ibx > 0) { errx = 0.06761; xrec -= 0.03832; }
3391 else { errx = 0.06653; xrec += 0.02581; }
3394 /* ---> Overflows */
3396 if (TMath::Abs(wx - 6) < 0.1) errx = 0.06979;
3397 else if (TMath::Abs(wx - 7.5) < 0.1) errx = 0.1089;
3398 else if (TMath::Abs(wx - 10) < 0.1) errx = 0.09847;
3402 if (TMath::Abs(wx - 6) < 0.1) errx = 0.06022;
3403 else if (TMath::Abs(wx - 7.5) < 0.1) errx = 0.07247;
3404 else if (TMath::Abs(wx - 10) < 0.1) errx = 0.07359;