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 Revision 1.16 2001/10/23 13:03:35 hristov
19 The access to several data members was changed from public to protected. The digitisation was adapted to the multi-event case (J.Chudoba)
21 Revision 1.15 2001/10/21 18:31:23 hristov
22 Several pointers were set to zero in the default constructors to avoid memory management problems
24 Revision 1.14 2001/05/14 13:25:54 hristov
25 stdlib.h included (for Alpha)
27 Revision 1.13 2001/05/10 12:26:31 jbarbosa
28 Totally reworked version of reconstruction algorithm.
30 Revision 1.12 2001/02/27 22:15:03 jbarbosa
31 Removed compiler warning.
33 Revision 1.11 2001/02/27 15:21:46 jbarbosa
34 Transition to SDigits.
36 Revision 1.10 2001/02/13 20:39:06 jbarbosa
37 Changes to make it work with new IO.
39 Revision 1.9 2001/01/22 21:39:11 jbarbosa
42 Revision 1.8 2000/11/15 15:52:53 jbarbosa
43 Turned on spot algorithm.
45 Revision 1.7 2000/11/01 15:37:05 jbarbosa
46 Updated to use its own rec. point object.
48 Revision 1.6 2000/10/02 21:28:12 fca
49 Removal of useless dependecies via forward declarations
51 Revision 1.5 2000/06/30 16:30:28 dibari
52 Disabled writing to rechits.
54 Revision 1.4 2000/06/15 15:46:59 jbarbosa
55 Corrected compilation errors on HP-UX (replaced pow with TMath::Power)
57 Revision 1.3 2000/06/13 13:15:41 jbarbosa
58 Still some code cleanup done (variable names)
60 Revision 1.2 2000/06/12 15:19:30 jbarbosa
63 Revision 1.1 2000/04/19 13:05:14 morsch
64 J. Barbosa's spot reconstruction algorithm.
72 #include "AliRICHPoints.h"
73 #include "AliRICHDetect.h"
74 #include "AliRICHHit.h"
75 #include "AliRICHDigit.h"
76 #include "AliRICHRawCluster.h"
77 #include "AliRICHCerenkov.h"
78 #include "AliRICHSegmentationV0.h"
80 #include "TParticle.h"
93 ClassImp(AliRICHDetect)
94 //___________________________________________
95 AliRICHDetect::AliRICHDetect() : TObject()
98 // Default constructor
106 //___________________________________________
107 AliRICHDetect::AliRICHDetect(const char *name, const char *title)
111 TStyle *mystyle=new TStyle("Plain","mystyle");
112 mystyle->SetPalette(1,0);
116 fc1= new TCanvas("c1","Reconstructed points",50,50,300,350);
118 fc2= new TCanvas("c2","Reconstructed points after SPOT",370,50,300,350);
120 fc3= new TCanvas("c3","Used Digits",690,50,300,350);
121 fc4= new TCanvas("c4","Mesh activation data",50,430,600,350);
127 //___________________________________________
128 AliRICHDetect::~AliRICHDetect()
136 void AliRICHDetect::Detect(Int_t nev, Int_t type)
140 // Detection algorithm
143 //printf("Detection started!\n");
144 Float_t omega,omega1,theta1,phi_relative,steptheta,stepphi,x,y,q=0,z,cx,cy,l,aux1,aux2,aux3,max,radius=0,meanradius=0;
145 Int_t maxi,maxj,maxk;
146 Float_t originalOmega, originalPhi, originalTheta;
147 Float_t binomega, bintheta, binphi;
148 Int_t intomega, inttheta, intphi;
151 AliRICH *pRICH = (AliRICH*)gAlice->GetDetector("RICH");
152 AliRICHSegmentationV0* segmentation;
153 AliRICHChamber* iChamber;
154 AliRICHGeometry* geometry;
156 iChamber = &(pRICH->Chamber(0));
157 segmentation=(AliRICHSegmentationV0*) iChamber->GetSegmentationModel(0);
158 geometry=iChamber->GetGeometryModel();
161 //const Float_t Noise_Level=0; //Noise Level in percentage of mesh points
162 //const Float_t t=0.6; //Softening of Noise Correction (factor)
164 const Float_t kPi=TMath::Pi();
166 const Float_t kHeight=geometry->GetRadiatorToPads(); //Distance from Radiator to Pads in centimeters
167 //printf("Distance to Pads:%f\n",kHeight);
169 const Int_t kSpot=0; //number of passes with spot algorithm
171 const Int_t kDimensionTheta=100; //Matrix dimension for angle Detection
172 const Int_t kDimensionPhi=100;
173 const Int_t kDimensionOmega=100;
175 const Float_t SPOTp=.25; //Percentage of spot action
176 const Float_t kMinOmega=.4;
177 const Float_t kMaxOmega=.7; //Maximum Cherenkov angle to identify
178 const Float_t kMinTheta=0;
179 const Float_t kMaxTheta=10*kPi/180;
180 const Float_t kMinPhi=0;
181 const Float_t kMaxPhi=360*kPi/180;
183 Float_t rechit[6]; //Reconstructed point data
185 Int_t ***point = i3tensor(0,kDimensionTheta,0,kDimensionPhi,0,kDimensionOmega);
186 Int_t ***point1 = i3tensor(0,kDimensionTheta,0,kDimensionPhi,0,kDimensionOmega);
188 steptheta=(kMaxTheta-kMinTheta)/kDimensionTheta;
189 stepphi=(kMaxPhi-kMinPhi)/kDimensionPhi;
191 static TH3F *Points = new TH3F("Points","Reconstructed points 3D",kDimensionTheta,0,kDimensionTheta,kDimensionPhi,0,kDimensionPhi,kDimensionOmega,0,kDimensionOmega);
192 static TH2F *ThetaPhi = new TH2F("ThetaPhi","Theta-Phi projection",kDimensionTheta,0,kDimensionTheta,kDimensionPhi,0,kDimensionPhi);
193 static TH2F *OmegaTheta = new TH2F("OmegaTheta","Omega-Theta projection",kDimensionTheta,0,kDimensionTheta,kDimensionOmega,0,kDimensionOmega);
194 static TH2F *OmegaPhi = new TH2F("OmegaPhi","Omega-Phi projection",kDimensionPhi,0,kDimensionPhi,kDimensionOmega,0,kDimensionOmega);
195 static TH3F *SpotPoints = new TH3F("Points","Reconstructed points 3D, spot",kDimensionTheta,0,kDimensionTheta,kDimensionPhi,0,kDimensionPhi,kDimensionOmega,0,kDimensionOmega);
196 static TH2F *SpotThetaPhi = new TH2F("ThetaPhi","Theta-Phi projection, spot",kDimensionTheta,0,kDimensionTheta,kDimensionPhi,0,kDimensionPhi);
197 static TH2F *SpotOmegaTheta = new TH2F("OmegaTheta","Omega-Theta projection, spot",kDimensionTheta,0,kDimensionTheta,kDimensionOmega,0,kDimensionOmega);
198 static TH2F *SpotOmegaPhi = new TH2F("OmegaPhi","Omega-Phi projection, spot",kDimensionPhi,0,kDimensionPhi,kDimensionOmega,0,kDimensionOmega);
199 static TH2F *DigitsXY = new TH2F("DigitsXY","Pads used for reconstruction",150,-25,25,150,-25,25);
200 static TH1F *AngleAct = new TH1F("AngleAct","Activation per angle",100,.45,1);
201 static TH1F *Activation = new TH1F("Activation","Activation per ring",100,0,25);
202 Points->SetXTitle("theta");
203 Points->SetYTitle("phi");
204 Points->SetZTitle("omega");
205 ThetaPhi->SetXTitle("theta");
206 ThetaPhi->SetYTitle("phi");
207 OmegaTheta->SetXTitle("theta");
208 OmegaTheta->SetYTitle("omega");
209 OmegaPhi->SetXTitle("phi");
210 OmegaPhi->SetYTitle("omega");
211 SpotPoints->SetXTitle("theta");
212 SpotPoints->SetYTitle("phi");
213 SpotPoints->SetZTitle("omega");
214 SpotThetaPhi->SetXTitle("theta");
215 SpotThetaPhi->SetYTitle("phi");
216 SpotOmegaTheta->SetXTitle("theta");
217 SpotOmegaTheta->SetYTitle("omega");
218 SpotOmegaPhi->SetXTitle("phi");
219 SpotOmegaPhi->SetYTitle("omega");
220 AngleAct->SetFillColor(5);
221 AngleAct->SetXTitle("rad");
222 AngleAct->SetYTitle("activation");
223 Activation->SetFillColor(5);
224 Activation->SetXTitle("activation");
226 Int_t ntracks = (Int_t)gAlice->TreeH()->GetEntries();
228 Float_t trackglob[3];
231 //printf("Area de uma elipse com teta 0 e Omega 45:%f",Area(0,45));
235 for (track=0; track<ntracks;track++) {
237 gAlice->TreeH()->GetEvent(track);
238 TClonesArray *pHits = pRICH->Hits();
239 if (pHits == 0) return;
240 Int_t nhits = pHits->GetEntriesFast();
241 if (nhits == 0) continue;
242 //Int_t nent=(Int_t)gAlice->TreeD()->GetEntries();
243 AliRICHHit *mHit = 0;
246 Int_t counter=0, counter1=0;
248 for(i=0;i<kDimensionTheta;i++)
250 for(j=0;j<kDimensionPhi;j++)
252 for(k=0;k<kDimensionOmega;k++)
256 //printf("Dimensions theta:%d, phi:%d, omega:%d",kDimensionTheta,kDimensionPhi,kDimensionOmega);
257 //printf("Resetting %d %d %d, time %d\n",i,j,k,counter);
258 //-Noise_Level*(Area(i*kPi/(18*dimension),k*kMaxOmega/dimension)-Area((i-1)*kPi/(18*dimension),(k-1)*kMaxOmega/dimension));
259 //printf("n-%f",-Noise_Level*(Area(i*kPi/(18*dimension),k*kMaxOmega/dimension)-Area((i-1)*kPi/(18*dimension),(k-1)*kMaxOmega/dimension)));
264 Int_t ncerenkovs = pRICH->Cerenkovs()->GetEntriesFast();
270 for (Int_t hit=0;hit<ncerenkovs;hit++) {
271 AliRICHCerenkov* cHit = (AliRICHCerenkov*) pRICH->Cerenkovs()->UncheckedAt(hit);
272 Float_t loss = cHit->fLoss; //did it hit the CsI?
273 Float_t production = cHit->fProduction; //was it produced in freon?
274 Float_t cherenkov = cHit->fCerenkovAngle; //production cerenkov angle
275 if (loss == 4 && production == 1)
278 originalOmega += cherenkov;
279 //printf("%f\n",cherenkov);
283 printf("Cerenkovs : %d\n",counter);
285 if(counter != 0) //if there are cerenkovs
287 originalOmega = originalOmega/counter;
288 printf("Original omega: %f\n",originalOmega);
292 mHit = (AliRICHHit*) pHits->UncheckedAt(0);
293 Int_t nch = mHit->Chamber();
294 originalTheta = mHit->Theta();
295 originalPhi = mHit->Phi();
296 trackglob[0] = mHit->X();
297 trackglob[1] = mHit->Y();
298 trackglob[2] = mHit->Z();
300 printf("\n--------------------------------------\n");
301 printf("Chamber %d, track %d\n", nch, track);
304 iChamber = &(pRICH->Chamber(nch-1));
306 //printf("Nch:%d\n",nch);
308 iChamber->GlobaltoLocal(trackglob,trackloc);
310 iChamber->LocaltoGlobal(trackloc,trackglob);
316 AliRICHDigit *points = 0;
317 TClonesArray *pDigits = pRICH->DigitsAddress(nch-1);
319 AliRICHRawCluster *cluster =0;
320 TClonesArray *pClusters = pRICH->RawClustAddress(nch-1);
324 //digitize from digits
328 gAlice->TreeD()->GetEvent(0);
329 Int_t ndigits = pDigits->GetEntriesFast();
331 printf("Got %d digits\n",ndigits);
334 //digitize from clusters
338 Int_t nent=(Int_t)gAlice->TreeR()->GetEntries();
339 gAlice->TreeR()->GetEvent(nent-1);
340 Int_t nclusters = pClusters->GetEntriesFast();
342 printf("Got %d clusters\n",nclusters);
349 printf("Starting calculations\n");
350 printf(" Start Finish\n");
351 printf("Progress: ");
352 for(Float_t theta=0;theta<kMaxTheta;theta+=steptheta)
355 for(Float_t phi=0;phi<=kMaxPhi;phi+=stepphi)
357 //printf("Phi:%3.1f\n", phi*180/kPi);
359 for (Int_t cycle=0;cycle<maxcycle;cycle++)
364 points=(AliRICHDigit*) pDigits->UncheckedAt(cycle);
365 segmentation->GetPadC(points->PadX(), points->PadY(),x, y, z);
370 cluster=(AliRICHRawCluster*) pClusters->UncheckedAt(cycle);
376 if(type ==0 || q > 100)
382 radius=TMath::Sqrt(TMath::Power(x,2)+TMath::Power(y,2));
384 //calculation of relative phi angle of digit
386 phi_relative = acos(y/radius);
387 phi_relative = TMath::Abs(phi_relative - phi);
396 //if (radius < (2*kHeight*tan(theta+kMaxOmega)))
397 if (radius < (2*kHeight*tan(kMaxOmega)))
398 //if(Fiducial(x,y,theta,phi,kHeight,kMaxOmega,kMinOmega))
403 //printf("Radius: %f, Max Radius: %f\n",radius,2*kHeight*tan(theta+kMaxOmega)*3/4);
404 //printf("Loaded digit %d with coordinates x:%f, y%f\n",dig,x,y);
405 //printf("Using digit %d, for theta %f\n",dig,theta);
410 l=kHeight/cos(theta);
414 DigitsXY->Fill(x,y,(float) 1);
416 theta1=SnellAngle(theta);
418 aux1=-y*sin(phi)+x*cos(phi);
419 aux2=y*cos(phi)+x*sin(phi);
420 aux3=( TMath::Power(aux1,2)+TMath::Power(cos(theta1)*aux2 ,2))/TMath::Power(sin(theta1)*aux2+l,2);
421 omega=atan(sqrt(aux3));
423 //omega is distorted, theta1 is distorted
425 if(InvSnellAngle(omega+TMath::Abs(cos(phi_relative))*theta1)<999)
427 omega1=InvSnellAngle(omega+TMath::Abs(cos(phi_relative))*theta);
428 theta1=InvSnellAngle(theta1);
437 if(omega1<kMaxOmega && omega1>kMinOmega)
439 //printf("Omega found:%f\n",omega);
440 omega1=omega1-kMinOmega;
442 //printf("Omega: %f Theta: %3.1f Phi:%3.1f\n",omega, theta*180/kPi, phi*180/kPi);
444 bintheta=theta1*kDimensionTheta/kMaxTheta;
445 binphi=phi*kDimensionPhi/kMaxPhi;
446 binomega=omega1*kDimensionOmega/(kMaxOmega-kMinOmega);
448 if(Int_t(bintheta+0.5)==Int_t(bintheta))
449 inttheta=Int_t(bintheta);
451 inttheta=Int_t(bintheta+0.5);
453 if(Int_t(binomega+0.5)==Int_t(binomega))
454 intomega=Int_t(binomega);
456 intomega=Int_t(binomega+0.5);
458 if(Int_t(binphi+0.5)==Int_t(binphi))
459 intphi=Int_t(binphi);
461 intphi=Int_t(binphi+0.5);
463 //printf("Point added at %d %d %d\n",inttheta,intphi,intomega);
466 point[inttheta][intphi][intomega]+=1;
468 point[inttheta][intphi][intomega]+=(Int_t)(q);
470 //printf("Omega stored:%d\n",intomega);
471 Points->Fill(inttheta,intphi,intomega,(float) 1);
472 ThetaPhi->Fill(inttheta,intphi,(float) 1);
473 OmegaTheta->Fill(inttheta,intomega,(float) 1);
474 OmegaPhi->Fill(intphi,intomega,(float) 1);
475 //printf("Filling at %d %d %d\n",Int_t(theta*kDimensionTheta/kMaxTheta),Int_t(phi*kDimensionPhi/kMaxPhi),Int_t(omega*kDimensionOmega/kMaxOmega));
477 //if(omega<kMaxOmega)point[Int_t(theta)][Int_t(phi)][Int_t(omega)]+=1;
483 //printf("Used %d digits for theta %3.1f\n",counter1, theta*180/kPi);
490 meanradius=meanradius/counter;
491 //printf("Mean radius:%f, counter:%d\n",meanradius,counter);
492 rechit[5]=meanradius;
493 //printf("Used %d digits\n",counter1);
501 Points->Draw("colz");
503 ThetaPhi->Draw("colz");
505 OmegaTheta->Draw("colz");
507 OmegaPhi->Draw("colz");
509 DigitsXY->Draw("colz");
514 //Points->Draw("same");
516 //ThetaPhi->Draw("same");
518 //OmegaTheta->Draw("same");
520 //OmegaPhi->Draw("same");
526 for(Int_t s=0;s<kSpot;s++)
528 printf(" Applying Spot algorithm, pass %d\n", s);
531 for(i=0;i<=kDimensionTheta;i++)
533 for(j=0;j<=kDimensionPhi;j++)
535 for(k=0;k<=kDimensionOmega;k++)
537 point1[i][j][k]=point[i][j][k];
543 for(i=1;i<kDimensionTheta-1;i++)
545 for(j=1;j<kDimensionPhi-1;j++)
547 for(k=1;k<kDimensionOmega-1;k++)
549 if((point[i][k][j]>point[i-1][k][j])&&(point[i][k][j]>point[i+1][k][j])&&
550 (point[i][k][j]>point[i][k-1][j])&&(point[i][k][j]>point[i][k+1][j])&&
551 (point[i][k][j]>point[i][k][j-1])&&(point[i][k][j]>point[i][k][j+1]))
553 //cout<<"SPOT"<<endl;
554 //Execute SPOT on point
555 point1[i][j][k]+=Int_t(SPOTp*(point[i-1][k][j]+point[i+1][k][j]+point[i][k-1][j]+point[i][k+1][j]+point[i][k][j-1]+point[i][k][j+1]));
556 point1[i-1][k][j]=Int_t(SPOTp*point[i-1][k][j]);
557 point1[i+1][k][j]=Int_t(SPOTp*point[i+1][k][j]);
558 point1[i][k-1][j]=Int_t(SPOTp*point[i][k-1][j]);
559 point1[i][k+1][j]=Int_t(SPOTp*point[i][k+1][j]);
560 point1[i][k][j-1]=Int_t(SPOTp*point[i][k][j-1]);
561 point1[i][k][j+1]=Int_t(SPOTp*point[i][k][j+1]);
567 //copy from buffer copy
569 for(i=1;i<kDimensionTheta;i++)
571 for(j=1;j<kDimensionPhi;j++)
573 for(k=1;k<kDimensionOmega;k++)
575 point[i][j][k]=point1[i][j][k];
580 if(point1[i][j][k] != 0)
582 SpotPoints->Fill(i,j,k,(float) point1[i][j][k]);
583 //printf("Random number %f\n",random->Rndm(2));
584 //if(random->Rndm() < .2)
586 SpotThetaPhi->Fill(i,j,(float) point1[i][j][k]);
587 SpotOmegaTheta->Fill(i,k,(float) point1[i][j][k]);
588 SpotOmegaPhi->Fill(j,k,(float) point1[i][j][k]);
591 //printf("Filling at %d %d %d value %f\n",i,j,k,(float) point1[i][j][k]);
595 //if(point1[i][j][k] != 0)
596 //printf("Last transfer point: %d, point1, %d\n",point[i][j][k],point1[i][j][k]);
602 //printf("Filled %d cells\n",counter1);
609 SpotPoints->Draw("colz");
611 SpotThetaPhi->Draw("colz");
613 SpotOmegaTheta->Draw("colz");
615 SpotOmegaPhi->Draw("colz");
620 //SpotPoints->Draw("same");
622 //SpotThetaPhi->Draw("same");
624 //SpotOmegaTheta->Draw("same");
626 //SpotOmegaPhi->Draw("same");
631 //Identification is equivalent to maximum determination
632 max=0;maxi=0;maxj=0;maxk=0;
634 printf(" Proceeding to identification");
636 for(i=0;i<kDimensionTheta;i++)
637 for(j=0;j<kDimensionPhi;j++)
638 for(k=0;k<kDimensionOmega;k++)
639 if(point[i][j][k]>max)
641 //cout<<"maxi="<<i*90/dimension<<" maxj="<<j*90/dimension<<" maxk="<<k*kMaxOmega/dimension*180/kPi<<" max="<<max<<endl;
642 maxi=i;maxj=j;maxk=k;
645 //printf("Max Omega %d, Max Theta %d, Max Phi %d (%d counts)\n",maxk,maxi,maxj,max);
649 Float_t FinalOmega = maxk*(kMaxOmega-kMinOmega)/kDimensionOmega;
650 Float_t FinalTheta = maxi*kMaxTheta/kDimensionTheta;
651 Float_t FinalPhi = maxj*kMaxPhi/kDimensionPhi;
653 FinalOmega += kMinOmega;
655 //printf("Detected angle for height %3.1f and for center %3.1f %3.1f:%f\n",h,cx,cy,maxk*kPi/(kDimensionTheta*4));
656 printf(" Indentified angles: cerenkov - %f, theta - %3.1f, phi - %3.1f (%f activation)\n", FinalOmega, FinalTheta*180/kPi, FinalPhi*180/kPi, max);
657 //printf("Detected angle for height %3.1f and for center %3.1f %3.1f:%f\n",kHeight,cx,cy,maxk);
659 AngleAct->Fill(FinalOmega, (float) max);
660 Activation->Fill(max, (float) 1);
672 AngleAct->Draw("same");
674 Activation->Draw("same");
678 //fscanf(omegas,"%f",&realomega);
679 //fscanf(thetas,"%f",&realtheta);
680 //printf("Real Omega: %f",realomega);
681 //cout<<"Detected:theta="<<maxi*90/kDimensionTheta<<"phi="<<maxj*90/kDimensionPhi<<"omega="<<maxk*kMaxOmega/kDimensionOmega*180/kPi<<" OmegaError="<<fabs(maxk*kMaxOmega/kDimensionOmega*180/kPi-realomega)<<" ThetaError="<<fabs(maxi*90/kDimensionTheta-realtheta)<<endl<<endl;
683 //fprintf(results,"Center Coordinates, cx=%6.2f cy=%6.2f, Real Omega=%6.2f, Detected Omega=%6.2f, Omega Error=%6.2f Theta Error=%6.2f\n",cx,cy,realomega,maxk*kMaxOmega/kDimensionOmega*180/kPi,fabs(maxk*kMaxOmega/kDimensionOmega*180/kPi-realomega),fabs(maxi*90/kDimensionTheta-realtheta));
686 pointpp(maxj*90/kDimensionTheta,maxi*90/kDimensionPhi,maxk*kMaxOmega/kDimensionOmega*180/kPi,cx,cy);//Generates a point on the elipse*/
689 //Start filling rec. hits
691 rechit[0] = FinalTheta;
692 rechit[1] = FinalPhi - 90*kPi/180;
693 rechit[2] = FinalOmega;
697 //CreatePoints(FinalTheta, 270*kPi/180 + FinalPhi, FinalOmega, kHeight);
699 printf ("track %d, theta %f, phi %f, omega %f\n\n\n",track,rechit[0],rechit[1]*180/kPi,rechit[2]);
702 pRICH->AddRecHit3D(nch-1,rechit,originalOmega, originalTheta, originalPhi);
703 //printf("rechit %f %f %f %f %f\n",rechit[0],rechit[1],rechit[2],rechit[3],rechit[4]);
704 //printf("Chamber:%d",nch);
707 else //if no cerenkovs
721 if(type==1) //reco from clusters
723 pRICH->ResetRawClusters();
724 //Int_t nent=(Int_t)gAlice->TreeR()->GetEntries();
725 //gAlice->TreeR()->GetEvent(track);
726 //printf("Going to branch %d\n",track);
727 //gAlice->GetEvent(nev);
731 //printf("\n\n\n\n");
732 gAlice->TreeR()->Fill();
734 for (i=0;i<kNCH;i++) {
735 fRec=pRICH->RecHitsAddress3D(i);
736 int ndig=fRec->GetEntriesFast();
737 printf ("Chamber %d, rings %d\n",i+1,ndig);
739 pRICH->ResetRecHits3D();
741 free_i3tensor(point,0,kDimensionTheta,0,kDimensionPhi,0,kDimensionOmega);
742 free_i3tensor(point1,0,kDimensionTheta,0,kDimensionPhi,0,kDimensionOmega);
747 Float_t AliRICHDetect:: Area(Float_t theta,Float_t omega)
751 // Calculates area of an ellipse for given incidence angles
755 const Float_t kHeight=9.25; //Distance from Radiator to Pads in pads
757 area=TMath::Pi()*TMath::Power(kHeight*tan(omega),2)/TMath::Power(TMath::Power(cos(theta),2)-TMath::Power(tan(omega)*sin(theta),2),3/2);
763 Int_t ***AliRICHDetect::i3tensor(long nrl, long nrh, long ncl, long nch, long ndl, long ndh)
764 // allocate a Int_t 3tensor with range t[nrl..nrh][ncl..nch][ndl..ndh]
766 long i,j,nrow=nrh-nrl+1,ncol=nch-ncl+1,ndep=ndh-ndl+1;
771 // allocate pointers to pointers to rows
772 t=(Int_t ***) malloc((size_t)((nrow+NR_END)*sizeof(Int_t**)));
773 if (!t) printf("allocation failure 1 in f3tensor()");
777 // allocate pointers to rows and set pointers to them
778 t[nrl]=(Int_t **) malloc((size_t)((nrow*ncol+NR_END)*sizeof(Int_t*)));
779 if (!t[nrl]) printf("allocation failure 2 in f3tensor()");
783 // allocate rows and set pointers to them
784 t[nrl][ncl]=(Int_t *) malloc((size_t)((nrow*ncol*ndep+NR_END)*sizeof(Int_t)));
785 if (!t[nrl][ncl]) printf("allocation failure 3 in f3tensor()");
786 t[nrl][ncl] += NR_END;
789 for(j=ncl+1;j<=nch;j++) t[nrl][j]=t[nrl][j-1]+ndep;
790 for(i=nrl+1;i<=nrh;i++) {
792 t[i][ncl]=t[i-1][ncl]+ncol*ndep;
793 for(j=ncl+1;j<=nch;j++) t[i][j]=t[i][j-1]+ndep;
796 // return pointer to array of pointers to rows
800 void AliRICHDetect::free_i3tensor(int ***t, long nrl, long nrh, long ncl, long nch,long ndl, long ndh)
801 // free a Int_t f3tensor allocated by i3tensor()
805 free((char*) (t[nrl][ncl]+ndl-NR_END));
806 free((char*) (t[nrl]+ncl-NR_END));
807 free((char*) (t+nrl-NR_END));
811 Float_t AliRICHDetect:: SnellAngle(Float_t iangle)
814 // Compute the Snell angle
816 Float_t nfreon = 1.295;
817 Float_t nquartz = 1.585;
827 sinrangle = a1*a2*sin(iangle);
834 rangle = asin(sinrangle);
838 Float_t AliRICHDetect:: InvSnellAngle(Float_t rangle)
841 // Compute the inverse Snell angle
843 Float_t nfreon = 1.295;
844 Float_t nquartz = 1.585;
854 siniangle = sin(rangle)/(a1*a2);
855 iangle = asin(siniangle);
862 iangle = asin(siniangle);
868 //________________________________________________________________________________
869 void AliRICHDetect::CreatePoints(Float_t theta, Float_t phi, Float_t omega, Float_t h)
872 // Create points along the ellipse equation
875 Float_t fiducial=h*TMath::Tan(omega+theta), l=h/TMath::Cos(theta), xtrial, y=0, c0, c1, c2;
876 //TRandom *random=new TRandom();
878 static TH2F *REllipse = new TH2F("REllipse","Reconstructed ellipses",150,-25,25,150,-25,25);
880 for(Float_t i=0;i<1000;i++)
886 while((c1*c1-4*c2*c0)<=0 && counter<1000)
888 //Choose which side to go...
889 if(i>250 && i<750) s1=1;
890 //if (gRandom->Rndm(1)>.5) s1=1;
892 //printf("s1:%d\n",s1);
894 y=s1*i*gRandom->Rndm(Int_t(fiducial/50));
895 //printf("Fiducial %f for omega:%f theta:%f phi:%f\n",fiducial,omega,theta,fphi);
898 Float_t omega1=omega;
900 //Solve the eq for a trial x
901 c0=-TMath::Power(y*TMath::Cos(alfa1)*TMath::Cos(theta1),2)-TMath::Power(y*TMath::Sin(alfa1),2)+TMath::Power(l*TMath::Tan(omega1),2)+2*l*y*TMath::Cos(alfa1)*TMath::Sin(theta1)*TMath::Power(TMath::Tan(omega1),2)+TMath::Power(y*TMath::Cos(alfa1)*TMath::Sin(theta1)*TMath::Tan(omega1),2);
902 c1=2*y*TMath::Cos(alfa1)*TMath::Sin(alfa1)-2*y*TMath::Cos(alfa1)*TMath::Power(TMath::Cos(theta1),2)*TMath::Sin(alfa1)+2*l*TMath::Sin(alfa1)*TMath::Sin(theta1)*TMath::Power(TMath::Tan(omega1),2)+2*y*TMath::Cos(alfa1)*TMath::Sin(alfa1)*TMath::Power(TMath::Sin(theta1),2)*TMath::Power(TMath::Tan(omega1),2);
903 c2=-TMath::Power(TMath::Cos(alfa1),2)-TMath::Power(TMath::Cos(theta1)*TMath::Sin(alfa1),2)+TMath::Power(TMath::Sin(alfa1)*TMath::Sin(theta1)*TMath::Tan(omega1),2);
904 //cout<<"Trial: y="<<y<<"c0="<<c0<<" c1="<<c1<<" c2="<<c2<<endl;
905 //printf("Result:%f\n\n",c1*c1-4*c2*c0);
913 //Choose which side to go...
914 //if(gRandom->Rndm(1)>.5) s=1;
917 //if (gRandom->Rndm(1)>.5) s2=1;
919 xtrial=(-c1+s2*TMath::Sqrt(c1*c1-4*c2*c0))/(2*c2);
920 //cout<<"x="<<xtrial<<" y="<<cy+y<<endl;
921 //printf("Coordinates: %f %f\n",xtrial,fCy+y);
923 REllipse->Fill(xtrial,y);
925 //printf("Coordinates: %f %f %f\n",vectorGlob[0],vectorGlob[1],vectorGlob[2]);
932 Int_t AliRICHDetect::Fiducial(Float_t rotx, Float_t roty, Float_t theta, Float_t phi, Float_t height, Float_t maxOmega, Float_t minOmega)
935 Float_t x,y,y1,h,omega1,omega2;
938 Float_t a1,b1, offset1;
942 //refraction calculation
944 theta = SnellAngle(theta);
945 phi = phi - TMath::Pi();
946 omega1 = SnellAngle(maxOmega);
947 omega2 = SnellAngle(minOmega);
950 a = h*(tan(omega1+theta)+tan(omega1-theta))/2;
953 offset = h*(tan(omega1+theta)-tan(omega1-theta))/2;
957 a1 = h*(tan(omega2+theta)+tan(omega2-theta))/2;
960 offset1 = h*(tan(omega2+theta)-tan(omega2-theta))/2;
964 x = rotx*cos(phi)+roty*sin(phi);
965 y = -rotx*sin(phi)+roty*cos(phi) - offset;
966 y1 = -rotx*sin(phi)+roty*cos(phi) - offset1;
969 if(x*x/a+y*y/b<1 && x*x/a1+y1*y1/b1>1)