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 #include "AliRICHClusterFinder.h"
21 #include "AliRICHHitMapA1.h"
22 #include "AliRICHSDigit.h"
23 #include "AliRICHDigit.h"
24 #include "AliRICHRawCluster.h"
32 #include <TPostScript.h>
35 //----------------------------------------------------------
36 static AliSegmentation* gSegmentation;
37 static AliRICHResponse* gResponse;
38 static Int_t gix[500];
39 static Int_t giy[500];
40 static Float_t gCharge[500];
42 static Int_t gFirst=kTRUE;
43 static TMinuit *gMyMinuit ;
44 void fcn(Int_t &npar, Double_t *gin, Double_t &f, Double_t *par, Int_t iflag);
45 static Int_t gChargeTot;
47 ClassImp(AliRICHClusterFinder)
49 AliRICHClusterFinder::AliRICHClusterFinder
50 (AliSegmentation *segmentation, AliRICHResponse *response,
51 TClonesArray *digits, Int_t chamber)
54 // Constructor for Cluster Finder object
56 fSegmentation=segmentation;
60 fNdigits = fDigits->GetEntriesFast();
62 fRawClusters=new TClonesArray("AliRICHRawCluster",10000);
71 AliRICHClusterFinder::AliRICHClusterFinder()
74 // Default constructor
82 fRawClusters=new TClonesArray("AliRICHRawCluster",10000);
92 AliRICHClusterFinder::AliRICHClusterFinder(const AliRICHClusterFinder& ClusterFinder)
93 :TObject(ClusterFinder)
98 AliRICHClusterFinder::~AliRICHClusterFinder()
106 void AliRICHClusterFinder::AddRawCluster(const AliRICHRawCluster c)
109 // Add a raw cluster copy to the list
111 AliRICH *pRICH=(AliRICH*)gAlice->GetModule("RICH");
112 pRICH->AddRawCluster(fChamber,c);
118 void AliRICHClusterFinder::Decluster(AliRICHRawCluster *cluster)
122 // Decluster algorithm
124 Int_t mul = cluster->fMultiplicity;
125 // printf("Decluster - multiplicity %d \n",mul);
127 if (mul == 1 || mul ==2) {
129 // Nothing special for 1- and 2-clusters
131 cluster->fNcluster[0]=fNPeaks;
132 cluster->fNcluster[1]=0;
134 AddRawCluster(*cluster);
136 } else if (mul ==3) {
138 // 3-cluster, check topology
139 // printf("\n 3-cluster, check topology \n");
140 if (fDeclusterFlag) {
141 if (Centered(cluster)) {
142 // ok, cluster is centered
144 // cluster is not centered, split into 2+1
148 cluster->fNcluster[0]=fNPeaks;
149 cluster->fNcluster[1]=0;
151 AddRawCluster(*cluster);
156 // 4-and more-pad clusters
158 if (mul <= fClusterSize) {
159 if (fDeclusterFlag) {
160 SplitByLocalMaxima(cluster);
163 cluster->fNcluster[0]=fNPeaks;
164 cluster->fNcluster[1]=0;
166 AddRawCluster(*cluster);
174 Bool_t AliRICHClusterFinder::Centered(AliRICHRawCluster *cluster)
177 // Is the cluster centered?
180 dig= (AliRICHDigit*)fDigits->UncheckedAt(cluster->fIndexMap[0]);
181 Int_t ix=dig->PadX();
182 Int_t iy=dig->PadY();
184 Int_t x[kMaxNeighbours], y[kMaxNeighbours], xN[kMaxNeighbours], yN[kMaxNeighbours];
186 fSegmentation->Neighbours(ix,iy,&nn,x,y);
188 for (Int_t i=0; i<nn; i++) {
189 if (fHitMap->TestHit(x[i],y[i]) == kUsed) {
194 //printf("Getting: %d %d %d\n",i,x[i],y[i]);
199 // cluster is centered !
201 cluster->fNcluster[0]=fNPeaks;
202 cluster->fNcluster[1]=0;
205 AddRawCluster(*cluster);
210 // Highest signal on an edge, split cluster into 2+1
212 // who is the neighbour ?
214 //printf("Calling GetIndex with x:%d y:%d\n",xN[0], yN[0]);
216 Int_t nind=fHitMap->GetHitIndex(xN[0], yN[0]);
217 Int_t i1= (nind==cluster->fIndexMap[1]) ? 1:2;
218 Int_t i2= (nind==cluster->fIndexMap[1]) ? 2:1;
221 AliRICHRawCluster cnew;
223 cnew.fNcluster[0]=-1;
224 cnew.fNcluster[1]=fNRawClusters;
226 cnew.fNcluster[0]=fNPeaks;
229 cnew.fMultiplicity=2;
230 cnew.fIndexMap[0]=cluster->fIndexMap[0];
231 cnew.fIndexMap[1]=cluster->fIndexMap[i1];
233 cnew.fClusterType=cnew.PhysicsContribution();
238 cluster->fMultiplicity=1;
239 cluster->fIndexMap[0]=cluster->fIndexMap[i2];
240 cluster->fIndexMap[1]=0;
241 cluster->fIndexMap[2]=0;
242 FillCluster(cluster);
244 cluster->fNcluster[0]=fNPeaks;
245 cluster->fNcluster[1]=0;
247 cluster->fClusterType=cluster->PhysicsContribution();
248 AddRawCluster(*cluster);
252 printf("\n Completely screwed up %d !! \n",nd);
258 void AliRICHClusterFinder::SplitByLocalMaxima(AliRICHRawCluster *c)
262 // Split the cluster according to the number of maxima inside
265 AliRICHDigit* dig[100], *digt;
266 Int_t ix[100], iy[100], q[100];
267 Float_t x[100], y[100], zdum;
268 Int_t i; // loops over digits
269 Int_t j; // loops over local maxima
272 // Int_t threshold=500;
273 Int_t mul=c->fMultiplicity;
275 // dump digit information into arrays
277 for (i=0; i<mul; i++)
279 dig[i]= (AliRICHDigit*)fDigits->UncheckedAt(c->fIndexMap[i]);
280 ix[i]= dig[i]->PadX();
281 iy[i]= dig[i]->PadY();
282 q[i] = dig[i]->Signal();
283 fSegmentation->GetPadC(ix[i], iy[i], x[i], y[i], zdum);
290 Int_t associatePeak[100];
293 Int_t xNei[kMaxNeighbours], yNei[kMaxNeighbours];
294 for (i=0; i<mul; i++) {
295 fSegmentation->Neighbours(ix[i], iy[i], &nn, xNei, yNei);
297 for (j=0; j<nn; j++) {
298 if (fHitMap->TestHit(xNei[j], yNei[j])==kEmpty) continue;
299 digt=(AliRICHDigit*) fHitMap->GetHit(xNei[j], yNei[j]);
300 if (digt->Signal() > q[i]) {
304 // handle special case of neighbouring pads with equal signal
305 } else if (digt->Signal() == q[i]) {
307 for (Int_t k=0; k<nLocal; k++) {
308 if (xNei[j]==ix[indLocal[k]] && yNei[j]==iy[indLocal[k]]){
314 } // loop over next neighbours
315 // Maxima should not be on the edge
320 } // loop over all digits
321 // printf("Found %d local Maxima",nLocal);
323 // If only one local maximum found but multiplicity is high
324 // take global maximum from the list of digits.
325 if (nLocal==1 && mul>5) {
327 for (i=0; i<mul; i++) {
338 // If number of local maxima is 2 try to fit a double gaussian
341 // Initialise global variables for fit
343 gSegmentation=fSegmentation;
344 gResponse =fResponse;
347 for (i=0; i<mul; i++) {
350 gCharge[i]=Float_t(q[i]);
355 gMyMinuit = new TMinuit(5);
357 gMyMinuit->SetFCN(fcn);
358 gMyMinuit->mninit(5,10,7);
359 Double_t arglist[20];
362 // gMyMinuit->mnexcm("SET ERR",arglist,1,ierflag);
363 // Set starting values
364 static Double_t vstart[5];
365 vstart[0]=x[indLocal[0]];
366 vstart[1]=y[indLocal[0]];
367 vstart[2]=x[indLocal[1]];
368 vstart[3]=y[indLocal[1]];
369 vstart[4]=Float_t(q[indLocal[0]])/
370 Float_t(q[indLocal[0]]+q[indLocal[1]]);
371 // lower and upper limits
372 static Double_t lower[5], upper[5];
373 Int_t isec=fSegmentation->Sector(ix[indLocal[0]], iy[indLocal[0]]);
374 lower[0]=vstart[0]-fSegmentation->Dpx(isec)/2;
375 lower[1]=vstart[1]-fSegmentation->Dpy(isec)/2;
376 // lower[1]=vstart[1];
378 upper[0]=lower[0]+fSegmentation->Dpx(isec);
379 upper[1]=lower[1]+fSegmentation->Dpy(isec);
380 // upper[1]=vstart[1];
382 isec=fSegmentation->Sector(ix[indLocal[1]], iy[indLocal[1]]);
383 lower[2]=vstart[2]-fSegmentation->Dpx(isec)/2;
384 lower[3]=vstart[3]-fSegmentation->Dpy(isec)/2;
385 // lower[3]=vstart[3];
387 upper[2]=lower[2]+fSegmentation->Dpx(isec);
388 upper[3]=lower[3]+fSegmentation->Dpy(isec);
389 // upper[3]=vstart[3];
394 static Double_t step[5]={0.005, 0.03, 0.005, 0.03, 0.01};
396 gMyMinuit->mnparm(0,"x1",vstart[0],step[0],lower[0],upper[0],ierflag);
397 gMyMinuit->mnparm(1,"y1",vstart[1],step[1],lower[1],upper[1],ierflag);
398 gMyMinuit->mnparm(2,"x2",vstart[2],step[2],lower[2],upper[2],ierflag);
399 gMyMinuit->mnparm(3,"y2",vstart[3],step[3],lower[3],upper[3],ierflag);
400 gMyMinuit->mnparm(4,"a0",vstart[4],step[4],lower[4],upper[4],ierflag);
401 // ready for minimisation
402 gMyMinuit->SetPrintLevel(-1);
403 gMyMinuit->mnexcm("SET OUT", arglist, 0, ierflag);
407 gMyMinuit->mnexcm("SET NOGR", arglist, 0, ierflag);
408 gMyMinuit->mnexcm("SCAN", arglist, 0, ierflag);
409 gMyMinuit->mnexcm("EXIT" , arglist, 0, ierflag);
411 // Double_t amin,edm,errdef;
412 // Int_t nvpar,nparx,icstat;
413 // gMyMinuit->mnstat(amin,edm,errdef,nvpar,nparx,icstat);
414 // gMyMinuit->mnprin(3,amin);
415 // Get fitted parameters
417 Double_t xrec[2], yrec[2], qfrac;
419 Double_t epxz, b1, b2;
421 gMyMinuit->mnpout(0, chname, xrec[0], epxz, b1, b2, ierflg);
422 gMyMinuit->mnpout(1, chname, yrec[0], epxz, b1, b2, ierflg);
423 gMyMinuit->mnpout(2, chname, xrec[1], epxz, b1, b2, ierflg);
424 gMyMinuit->mnpout(3, chname, yrec[1], epxz, b1, b2, ierflg);
425 gMyMinuit->mnpout(4, chname, qfrac, epxz, b1, b2, ierflg);
426 //printf("\n %f %f %f %f %f\n", xrec[0], yrec[0], xrec[1], yrec[1],qfrac);
431 // One cluster for each maximum
433 for (j=0; j<2; j++) {
434 AliRICHRawCluster cnew;
436 cnew.fNcluster[0]=-1;
437 cnew.fNcluster[1]=fNRawClusters;
439 cnew.fNcluster[0]=fNPeaks;
442 cnew.fMultiplicity=0;
443 cnew.fX=Float_t(xrec[j]);
444 cnew.fY=Float_t(yrec[j]);
446 cnew.fQ=Int_t(gChargeTot*qfrac);
448 cnew.fQ=Int_t(gChargeTot*(1-qfrac));
450 gSegmentation->SetHit(xrec[j],yrec[j],0);
451 for (i=0; i<mul; i++) {
452 cnew.fIndexMap[cnew.fMultiplicity]=c->fIndexMap[i];
453 gSegmentation->SetPad(gix[i], giy[i]);
454 Float_t q1=gResponse->IntXY(gSegmentation);
455 cnew.fContMap[cnew.fMultiplicity]=Float_t(q[i])/(q1*cnew.fQ);
456 cnew.fMultiplicity++;
458 FillCluster(&cnew,0);
459 //printf("\n x,y %f %f ", cnew.fX, cnew.fY);
460 cnew.fClusterType=cnew.PhysicsContribution();
468 if (nLocal !=-100 || !fitted) {
469 // Check if enough local clusters have been found,
470 // if not add global maxima to the list
476 printf("\n Warning, no local maximum found \n");
480 if (nPerMax > fNperMax) {
481 Int_t nGlob=mul/fNperMax-nLocal+1;
484 for (i=0; i<mul; i++) {
491 if (nnew==nGlob) break;
496 // Associate hits to peaks
498 for (i=0; i<mul; i++) {
501 if (isLocal[i]) continue;
502 for (j=0; j<nLocal; j++) {
503 Int_t il=indLocal[j];
504 Float_t d=TMath::Sqrt((x[i]-x[il])*(x[i]-x[il])
505 +(y[i]-y[il])*(y[i]-y[il]));
508 // Select nearest peak
514 } else if (d==dmin) {
516 // If more than one take highest peak
529 // One cluster for each maximum
531 for (j=0; j<nLocal; j++) {
532 AliRICHRawCluster cnew;
534 cnew.fNcluster[0]=-1;
535 cnew.fNcluster[1]=fNRawClusters;
537 cnew.fNcluster[0]=fNPeaks;
540 cnew.fIndexMap[0]=c->fIndexMap[indLocal[j]];
541 cnew.fMultiplicity=1;
542 for (i=0; i<mul; i++) {
543 if (isLocal[i]) continue;
544 if (associatePeak[i]==j) {
545 cnew.fIndexMap[cnew.fMultiplicity]=c->fIndexMap[i];
546 cnew.fMultiplicity++;
550 cnew.fClusterType=cnew.PhysicsContribution();
558 void AliRICHClusterFinder::FillCluster(AliRICHRawCluster* c, Int_t flag)
561 // Completes cluster information starting from list of digits
577 for (Int_t i=0; i<c->fMultiplicity; i++)
579 dig= (AliRICHDigit*)fDigits->UncheckedAt(c->fIndexMap[i]);
580 ix=dig->PadX()+c->fOffsetMap[i];
582 Int_t q=dig->Signal();
583 if (dig->Physics() >= dig->Signal()) {
585 } else if (dig->Physics() == 0) {
587 } else c->fPhysicsMap[i]=1;
590 // peak signal and track list
592 if (q>c->fPeakSignal) {
595 c->fTracks[0]=dig->Track(0);
596 c->fTracks[1]=dig->Track(1);
597 c->fTracks[2]=dig->Track(2);
599 //c->fTracks[0]=dig->fTrack;
600 c->fTracks[0]=dig->Hit();
601 c->fTracks[1]=dig->Track(0);
602 c->fTracks[2]=dig->Track(1);
605 if (c->fContMap[i] > frac) {
609 c->fTracks[0]=dig->Track(0);
610 c->fTracks[1]=dig->Track(1);
611 c->fTracks[2]=dig->Track(2);
613 //c->fTracks[0]=dig->fTrack;
614 c->fTracks[0]=dig->Hit();
615 c->fTracks[1]=dig->Track(0);
616 c->fTracks[2]=dig->Track(1);
621 fSegmentation->GetPadC(ix, iy, x, y, z);
627 } // loop over digits
632 c->fX=fSegmentation->GetAnod(c->fX);
635 // apply correction to the coordinate along the anode wire
639 fSegmentation->GetPadI(x, y, 0, ix, iy);
640 fSegmentation->GetPadC(ix, iy, x, y, z);
641 Int_t isec=fSegmentation->Sector(ix,iy);
642 TF1* cogCorr = fSegmentation->CorrFunc(isec-1);
645 Float_t yOnPad=(c->fY-y)/fSegmentation->Dpy(isec);
646 c->fY=c->fY-cogCorr->Eval(yOnPad, 0, 0);
652 void AliRICHClusterFinder::FindCluster(Int_t i, Int_t j, AliRICHRawCluster &c){
657 // Add i,j as element of the cluster
660 Int_t idx = fHitMap->GetHitIndex(i,j);
661 AliRICHDigit* dig = (AliRICHDigit*) fHitMap->GetHit(i,j);
662 Int_t q=dig->Signal();
663 if (q > TMath::Abs(c.fPeakSignal)) {
666 c.fTracks[0]=dig->fTracks[0];
667 c.fTracks[1]=dig->fTracks[1];
668 c.fTracks[2]=dig->fTracks[2];
670 //c.fTracks[0]=dig->fTrack;
671 c.fTracks[0]=dig->Hit();
672 c.fTracks[1]=dig->Track(0);
673 c.fTracks[2]=dig->Track(1);
676 // Make sure that list of digits is ordered
678 Int_t mu=c.fMultiplicity;
681 if (dig->Physics() >= dig->Signal()) {
683 } else if (dig->Physics() == 0) {
685 } else c.fPhysicsMap[mu]=1;
688 for (Int_t ind=mu-1; ind>=0; ind--) {
689 Int_t ist=(c.fIndexMap)[ind];
690 Int_t ql=((AliRICHDigit*)fDigits
691 ->UncheckedAt(ist))->Signal();
693 c.fIndexMap[ind]=idx;
694 c.fIndexMap[ind+1]=ist;
703 if (c.fMultiplicity >= 50 ) {
704 printf("FindCluster - multiplicity >50 %d \n",c.fMultiplicity);
708 // Prepare center of gravity calculation
710 fSegmentation->GetPadC(i, j, x, y, z);
715 fHitMap->FlagHit(i,j);
717 // Now look recursively for all neighbours
720 Int_t xList[kMaxNeighbours], yList[kMaxNeighbours];
721 fSegmentation->Neighbours(i,j,&nn,xList,yList);
722 for (Int_t in=0; in<nn; in++) {
725 if (fHitMap->TestHit(ix,iy)==kUnused) FindCluster(ix, iy, c);
729 //_____________________________________________________________________________
731 void AliRICHClusterFinder::FindRawClusters()
734 // simple RICH cluster finder from digits -- finds neighbours and
735 // fill the tree with raw clusters
737 if (!fNdigits) return;
739 fHitMap = new AliRICHHitMapA1(fSegmentation, fDigits);
743 //printf ("Now I'm here");
749 for (ndig=0; ndig<fNdigits; ndig++) {
750 dig = (AliRICHDigit*)fDigits->UncheckedAt(ndig);
753 if (fHitMap->TestHit(i,j)==kUsed ||fHitMap->TestHit(i,j)==kEmpty) {
759 c.fPeakSignal=dig->Signal();
761 c.fTracks[0]=dig->fTracks[0];
762 c.fTracks[1]=dig->fTracks[1];
763 c.fTracks[2]=dig->fTracks[2];
765 //c.fTracks[0]=dig->fTrack;
766 c.fTracks[0]=dig->Hit();
767 c.fTracks[1]=dig->Track(0);
768 c.fTracks[2]=dig->Track(1);
769 // tag the beginning of cluster list in a raw cluster
774 c.fX=fSegmentation->GetAnod(c.fX);
777 // apply correction to the coordinate along the anode wire
784 fSegmentation->GetPadI(x, y, 0, ix, iy);
785 fSegmentation->GetPadC(ix, iy, x, y, z);
786 Int_t isec=fSegmentation->Sector(ix,iy);
787 TF1* cogCorr=fSegmentation->CorrFunc(isec-1);
789 Float_t yOnPad=(c.fY-y)/fSegmentation->Dpy(isec);
790 c.fY=c.fY-cogCorr->Eval(yOnPad,0,0);
794 // Analyse cluster and decluster if necessary
797 c.fNcluster[1]=fNRawClusters;
798 c.fClusterType=c.PhysicsContribution();
804 // reset Cluster object
805 for (int k=0;k<c.fMultiplicity;k++) {
813 void AliRICHClusterFinder::
824 fSegmentation->GiveTestPoints(n, x, y);
825 for (i=0; i<n; i++) {
829 SinoidalFit(xtest, ytest, func);
830 if (func) fSegmentation->SetCorrFunc(i, new TF1(*func));
836 void AliRICHClusterFinder::
837 SinoidalFit(Float_t x, Float_t y, TF1 *func)
842 static Int_t count=0;
847 sprintf(canvasname,"c%d",count);
850 Float_t xg[kNs], yg[kNs], xrg[kNs], yrg[kNs];
851 Float_t xsig[kNs], ysig[kNs];
853 AliSegmentation *segmentation=fSegmentation;
856 segmentation->GetPadI(x,y,0,ix,iy);
857 segmentation->GetPadC(ix,iy,x,y,z);
858 Int_t isec=segmentation->Sector(ix,iy);
860 Float_t xmin = x-segmentation->Dpx(isec)/2;
861 Float_t ymin = y-segmentation->Dpy(isec)/2;
863 // Integration Limits
864 Float_t dxI=fResponse->SigmaIntegration()*fResponse->ChargeSpreadX();
865 Float_t dyI=fResponse->SigmaIntegration()*fResponse->ChargeSpreadY();
875 Float_t dy=segmentation->Dpy(isec)/(kNs-1);
877 for (i=0; i<kNs; i++) {
883 segmentation->SigGenInit(x, yscan, 0);
885 for (segmentation->FirstPad(x, yscan,0, dxI, dyI);
886 segmentation->MorePads();
887 segmentation->NextPad())
889 qp=fResponse->IntXY(segmentation);
895 Int_t ixs=segmentation->Ix();
896 Int_t iys=segmentation->Iy();
898 segmentation->GetPadC(ixs,iys,xs,ys,zs);
902 Float_t ycog=sum/qcheck;
903 yg[i]=(yscan-y)/segmentation->Dpy(isec);
904 yrg[i]=(ycog-y)/segmentation->Dpy(isec);
911 Float_t dx=segmentation->Dpx(isec)/(kNs-1);
913 for (i=0; i<kNs; i++) {
919 segmentation->SigGenInit(xscan, y, 0);
921 for (segmentation->FirstPad(xscan, y, 0, dxI, dyI);
922 segmentation->MorePads();
923 segmentation->NextPad())
925 qp=fResponse->IntXY(segmentation);
931 Int_t ixs=segmentation->Ix();
932 Int_t iys=segmentation->Iy();
934 segmentation->GetPadC(ixs,iys,xs,ys,zs);
938 Float_t xcog=sum/qcheck;
939 xcog=segmentation->GetAnod(xcog);
941 xg[i]=(xscan-x)/segmentation->Dpx(isec);
942 xrg[i]=(xcog-x)/segmentation->Dpx(isec);
947 // Creates a Root function based on function sinoid above
948 // and perform the fit
950 // TGraph *graphx = new TGraph(kNs,xg ,xsig);
951 // TGraph *graphxr= new TGraph(kNs,xrg,xsig);
952 // TGraph *graphy = new TGraph(kNs,yg ,ysig);
953 TGraph *graphyr= new TGraph(kNs,yrg,ysig);
955 Double_t sinoid(Double_t *x, Double_t *par);
956 new TF1("sinoidf",sinoid,0.5,0.5,5);
957 graphyr->Fit("sinoidf","Q");
958 func = (TF1*)graphyr->GetListOfFunctions()->At(0);
961 TCanvas *c1=new TCanvas(canvasname,canvasname,400,10,600,700);
962 TPad* pad11 = new TPad("pad11"," ",0.01,0.51,0.49,0.99);
963 TPad* pad12 = new TPad("pad12"," ",0.51,0.51,0.99,0.99);
964 TPad* pad13 = new TPad("pad13"," ",0.01,0.01,0.49,0.49);
965 TPad* pad14 = new TPad("pad14"," ",0.51,0.01,0.99,0.49);
966 pad11->SetFillColor(11);
967 pad12->SetFillColor(11);
968 pad13->SetFillColor(11);
969 pad14->SetFillColor(11);
977 graphx->SetFillColor(42);
978 graphx->SetMarkerColor(4);
979 graphx->SetMarkerStyle(21);
981 graphx->GetHistogram()->SetXTitle("x on pad");
982 graphx->GetHistogram()->SetYTitle("xcog-x");
986 graphxr->SetFillColor(42);
987 graphxr->SetMarkerColor(4);
988 graphxr->SetMarkerStyle(21);
990 graphxr->GetHistogram()->SetXTitle("xcog on pad");
991 graphxr->GetHistogram()->SetYTitle("xcog-x");
995 graphy->SetFillColor(42);
996 graphy->SetMarkerColor(4);
997 graphy->SetMarkerStyle(21);
999 graphy->GetHistogram()->SetXTitle("y on pad");
1000 graphy->GetHistogram()->SetYTitle("ycog-y");
1005 graphyr->SetFillColor(42);
1006 graphyr->SetMarkerColor(4);
1007 graphyr->SetMarkerStyle(21);
1008 graphyr->Draw("AF");
1009 graphyr->GetHistogram()->SetXTitle("ycog on pad");
1010 graphyr->GetHistogram()->SetYTitle("ycog-y");
1016 Double_t sinoid(Double_t *x, Double_t *par)
1021 Double_t arg = -2*TMath::Pi()*x[0];
1022 Double_t fitval= par[0]*TMath::Sin(arg)+
1023 par[1]*TMath::Sin(2*arg)+
1024 par[2]*TMath::Sin(3*arg)+
1025 par[3]*TMath::Sin(4*arg)+
1026 par[4]*TMath::Sin(5*arg);
1031 Double_t DoubleGauss(Double_t *x, Double_t *par)
1034 // Doublr gaussian function
1036 Double_t arg1 = (x[0]-par[1])/0.18;
1037 Double_t arg2 = (x[0]-par[3])/0.18;
1038 Double_t fitval= par[0]*TMath::Exp(-arg1*arg1/2)
1039 +par[2]*TMath::Exp(-arg2*arg2/2);
1043 Float_t DiscrCharge(Int_t i,Double_t *par)
1045 // par[0] x-position of first cluster
1046 // par[1] y-position of first cluster
1047 // par[2] x-position of second cluster
1048 // par[3] y-position of second cluster
1049 // par[4] charge fraction of first cluster
1050 // 1-par[4] charge fraction of second cluster
1052 static Float_t qtot;
1055 for (Int_t jbin=0; jbin<gNbins; jbin++) {
1056 qtot+=gCharge[jbin];
1059 //printf("\n sum of charge from DiscrCharge %f\n", qtot);
1060 gChargeTot=Int_t(qtot);
1063 gSegmentation->SetPad(gix[i], giy[i]);
1065 gSegmentation->SetHit(par[0],par[1],0);
1066 Float_t q1=gResponse->IntXY(gSegmentation);
1069 gSegmentation->SetHit(par[2],par[3],0);
1070 Float_t q2=gResponse->IntXY(gSegmentation);
1072 Float_t value = qtot*(par[4]*q1+(1.-par[4])*q2);
1077 // Minimisation function
1078 void fcn(Int_t &npar, Double_t */*gin*/, Double_t &f, Double_t *par, Int_t /*iflag*/)
1087 for (i=0; i<gNbins; i++) {
1088 Float_t q0=gCharge[i];
1089 Float_t q1=DiscrCharge(i,par);
1090 delta=(q0-q1)/TMath::Sqrt(q0);
1095 chisq=chisq+=(qtot-qcont)*(qtot-qcont)*0.5;
1100 void AliRICHClusterFinder::SetDigits(TClonesArray *RICHdigits)
1103 // Get all the digits
1106 fNdigits = fDigits->GetEntriesFast();
1109 AliRICHClusterFinder& AliRICHClusterFinder::operator=(const AliRICHClusterFinder& /*rhs*/)
1111 // Assignment operator