]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RICH/AliRICHDetect.cxx
stdlib.h included (for Alpha)
[u/mrichter/AliRoot.git] / RICH / AliRICHDetect.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /*
17   $Log$
18   Revision 1.13  2001/05/10 12:26:31  jbarbosa
19   Totally reworked version of reconstruction algorithm.
20
21   Revision 1.12  2001/02/27 22:15:03  jbarbosa
22   Removed compiler warning.
23
24   Revision 1.11  2001/02/27 15:21:46  jbarbosa
25   Transition to SDigits.
26
27   Revision 1.10  2001/02/13 20:39:06  jbarbosa
28   Changes to make it work with new IO.
29
30   Revision 1.9  2001/01/22 21:39:11  jbarbosa
31   Several tune-ups
32
33   Revision 1.8  2000/11/15 15:52:53  jbarbosa
34   Turned on spot algorithm.
35
36   Revision 1.7  2000/11/01 15:37:05  jbarbosa
37   Updated to use its own rec. point object.
38
39   Revision 1.6  2000/10/02 21:28:12  fca
40   Removal of useless dependecies via forward declarations
41
42   Revision 1.5  2000/06/30 16:30:28  dibari
43   Disabled writing to rechits.
44
45   Revision 1.4  2000/06/15 15:46:59  jbarbosa
46   Corrected compilation errors on HP-UX (replaced pow with TMath::Power)
47
48   Revision 1.3  2000/06/13 13:15:41  jbarbosa
49   Still some code cleanup done (variable names)
50
51   Revision 1.2  2000/06/12 15:19:30  jbarbosa
52   Cleaned up version.
53
54   Revision 1.1  2000/04/19 13:05:14  morsch
55   J. Barbosa's spot reconstruction algorithm.
56
57 */
58
59 #include <stdlib.h>
60
61
62 #include "AliRICH.h"
63 #include "AliRICHPoints.h"
64 #include "AliRICHDetect.h"
65 #include "AliRICHHit.h"
66 #include "AliRICHDigit.h"
67 #include "AliRICHSegmentationV0.h"
68 #include "AliRun.h"
69 #include "TParticle.h"
70 #include "TTree.h"
71 #include "TMath.h"
72 #include "TRandom.h"
73 #include "TH3.h"
74 #include "TH2.h"
75 #include "TCanvas.h"
76
77 #include "malloc.h"
78
79
80 ClassImp(AliRICHDetect)
81 //___________________________________________
82 AliRICHDetect::AliRICHDetect() : TObject()
83 {
84
85 // Default constructor 
86
87 }
88
89 //___________________________________________
90 AliRICHDetect::AliRICHDetect(const char *name, const char *title)
91     : TObject()
92 {
93
94
95   fc1= new TCanvas("c1","Reconstructed points",50,50,300,350);
96   fc1->Divide(2,2);
97   fc2= new TCanvas("c2","Reconstructed points after SPOT",50,50,300,350);
98   fc2->Divide(2,2); 
99   fc3= new TCanvas("c3","Used Digits",50,50,300,350);
100   //fc3->Divide(2,1); 
101
102 }
103
104 //___________________________________________
105 AliRICHDetect::~AliRICHDetect()
106 {
107     
108 // Destructor
109
110 }
111
112
113 void AliRICHDetect::Detect(Int_t nev)
114 {       
115     
116 //
117 // Detection algorithm
118
119
120   //printf("Detection started!\n");
121   Float_t omega,omega1,theta1,steptheta,stepphi,x,y,z,cx,cy,l,aux1,aux2,aux3,max,radius=0,meanradius=0;
122   Int_t maxi,maxj,maxk;
123   //Float_t theta,phi,realomega,realtheta;
124   Float_t binomega, bintheta, binphi;
125   Int_t intomega, inttheta, intphi;
126   Int_t i,j,k;
127
128   AliRICH *pRICH  = (AliRICH*)gAlice->GetDetector("RICH");
129   AliRICHSegmentationV0*  segmentation;
130   AliRICHChamber*       iChamber;
131   AliRICHGeometry*  geometry;
132   
133   iChamber = &(pRICH->Chamber(0));
134   segmentation=(AliRICHSegmentationV0*) iChamber->GetSegmentationModel(0);
135   geometry=iChamber->GetGeometryModel();
136  
137   
138   //const Float_t Noise_Level=0;                       //Noise Level in percentage of mesh points
139   //const Float_t t=0.6;                               //Softening of Noise Correction (factor)
140   
141   const Float_t kPi=TMath::Pi();                
142   
143   const Float_t kHeight=geometry->GetRadiatorToPads(); //Distance from Radiator to Pads in centimeters
144   //printf("Distance to Pads:%f\n",kHeight);
145  
146   const Int_t kSpot=0;                                 //number of passes with spot algorithm
147   
148   const Int_t kDimensionTheta=30;                      //Matrix dimension for angle Detection
149   const Int_t kDimensionPhi=45;
150   const Int_t kDimensionOmega=100;
151   
152   const Float_t SPOTp=1;                              //Percentage of spot action
153   const Float_t kMinOmega=20*kPi/180;
154   const Float_t kMaxOmega=70*kPi/180;                 //Maximum Cherenkov angle to identify
155   const Float_t kMinTheta=0;
156   const Float_t kMaxTheta=15*kPi/180;   
157   //const Float_t kMaxTheta=0.1;
158   const Float_t kMinPhi=0;
159   const Float_t kMaxPhi=360*kPi/180;
160
161  
162   Float_t kCorr=0.61;                              //Correction factor, accounting for aberration, refractive index, etc.
163   //const Float_t kCorr=.9369;                        //from 0 incidence  
164   //const Float_t kCorr=1;
165
166   //TRandom* random=0;
167
168   Float_t rechit[6];                                 //Reconstructed point data
169
170   
171
172   //printf("Creating matrices\n");
173   //Float_t point[kDimensionTheta][kDimensionPhi][kDimensionOmega];
174   //Float_t point1[kDimensionTheta][kDimensionPhi][kDimensionOmega];
175   //printf("Created matrices\n");
176
177   Int_t ***point = i3tensor(0,kDimensionTheta,0,kDimensionPhi,0,kDimensionOmega);
178   Int_t ***point1 = i3tensor(0,kDimensionTheta,0,kDimensionPhi,0,kDimensionOmega);
179   
180   //Int_t **point  = new Int_t[kDimensionTheta][kDimensionPhi][kDimensionOmega];
181   //Int_t **point1 = new Int_t[kDimensionTheta][kDimensionPhi][kDimensionOmega];
182
183   steptheta=(kMaxTheta-kMinTheta)/kDimensionTheta;
184   stepphi=(kMaxPhi-kMinPhi)/kDimensionPhi;
185
186   static TH3F *Points = new TH3F("Points","Reconstructed points 3D",kDimensionTheta,0,kDimensionTheta,kDimensionPhi,0,kDimensionPhi,kDimensionOmega,0,kDimensionOmega);
187   static TH2F *ThetaPhi = new TH2F("ThetaPhi","Theta-Phi projection",kDimensionTheta,0,kDimensionTheta,kDimensionPhi,0,kDimensionPhi);
188   static TH2F *OmegaTheta = new TH2F("OmegaTheta","Omega-Theta projection",kDimensionTheta,0,kDimensionTheta,kDimensionOmega,0,kDimensionOmega);
189   static TH2F *OmegaPhi = new TH2F("OmegaPhi","Omega-Phi projection",kDimensionPhi,0,kDimensionPhi,kDimensionOmega,0,kDimensionOmega);
190   static TH3F *SpotPoints = new TH3F("Points","Reconstructed points 3D, spot",kDimensionTheta,0,kDimensionTheta,kDimensionPhi,0,kDimensionPhi,kDimensionOmega,0,kDimensionOmega);
191   static TH2F *SpotThetaPhi = new TH2F("ThetaPhi","Theta-Phi projection, spot",kDimensionTheta,0,kDimensionTheta,kDimensionPhi,0,kDimensionPhi);
192   static TH2F *SpotOmegaTheta = new TH2F("OmegaTheta","Omega-Theta projection, spot",kDimensionTheta,0,kDimensionTheta,kDimensionOmega,0,kDimensionOmega);
193   static TH2F *SpotOmegaPhi = new TH2F("OmegaPhi","Omega-Phi projection, spot",kDimensionPhi,0,kDimensionPhi,kDimensionOmega,0,kDimensionOmega);
194   static TH2F *DigitsXY = new TH2F("DigitsXY","Pads used for reconstruction",150,-25,25,150,-25,25);
195   Points->SetXTitle("theta");
196   Points->SetYTitle("phi");
197   Points->SetZTitle("omega");
198   ThetaPhi->SetXTitle("theta");
199   ThetaPhi->SetYTitle("phi");
200   OmegaTheta->SetXTitle("theta");
201   OmegaTheta->SetYTitle("omega");
202   OmegaPhi->SetXTitle("phi");
203   OmegaPhi->SetYTitle("omega");
204   SpotPoints->SetXTitle("theta");
205   SpotPoints->SetYTitle("phi");
206   SpotPoints->SetZTitle("omega");
207   SpotThetaPhi->SetXTitle("theta");
208   SpotThetaPhi->SetYTitle("phi");
209   SpotOmegaTheta->SetXTitle("theta");
210   SpotOmegaTheta->SetYTitle("omega");
211   SpotOmegaPhi->SetXTitle("phi");
212   SpotOmegaPhi->SetYTitle("omega");
213
214   Int_t ntracks = (Int_t)gAlice->TreeH()->GetEntries();
215   //Int_t ntrks = gAlice->GetNtrack();
216   
217   Float_t trackglob[3];
218   Float_t trackloc[3];
219
220   //printf("Area de uma elipse com teta 0 e Omega 45:%f",Area(0,45));
221     
222   Int_t track;
223         
224   for (track=0; track<ntracks;track++) {
225     gAlice->ResetHits();
226     gAlice->TreeH()->GetEvent(track);
227     TClonesArray *pHits  = pRICH->Hits();
228     if (pHits == 0) return;
229     Int_t nhits = pHits->GetEntriesFast();
230     if (nhits == 0) continue;
231     //Int_t nent=(Int_t)gAlice->TreeD()->GetEntries();
232     gAlice->TreeD()->GetEvent(0);
233     AliRICHHit *mHit = 0;
234     AliRICHDigit *points = 0;
235     //Int_t npoints=0;
236     
237     Int_t counter=0, counter1=0;
238     //Initialization
239     for(i=0;i<kDimensionTheta;i++)
240       {
241         for(j=0;j<kDimensionPhi;j++)
242           {
243             for(k=0;k<kDimensionOmega;k++)
244               {
245                 counter++;
246                 point[i][j][k]=0;
247                 //printf("Dimensions theta:%d, phi:%d, omega:%d",kDimensionTheta,kDimensionPhi,kDimensionOmega);
248                 //printf("Resetting %d %d %d, time %d\n",i,j,k,counter);
249                 //-Noise_Level*(Area(i*kPi/(18*dimension),k*kMaxOmega/dimension)-Area((i-1)*kPi/(18*dimension),(k-1)*kMaxOmega/dimension));
250                 //printf("n-%f",-Noise_Level*(Area(i*kPi/(18*dimension),k*kMaxOmega/dimension)-Area((i-1)*kPi/(18*dimension),(k-1)*kMaxOmega/dimension)));
251               }
252           }
253       }
254     mHit = (AliRICHHit*) pHits->UncheckedAt(0);
255     //printf("Aqui vou eu\n");
256     Int_t nch  = mHit->fChamber;
257     //printf("Aqui fui eu\n");
258     trackglob[0] = mHit->X();
259     trackglob[1] = mHit->Y();
260     trackglob[2] = mHit->Z();
261
262     printf("Chamber processed:%d\n",nch);
263
264     printf("Reconstructing particle at (global coordinates): %3.1f %3.1f %3.1f,\n",trackglob[0],trackglob[1],trackglob[2]);
265
266     iChamber = &(pRICH->Chamber(nch-1));
267     
268     //printf("Nch:%d\n",nch);
269
270     iChamber->GlobaltoLocal(trackglob,trackloc);
271     
272     printf("Reconstructing particle at (local coordinates) : %3.1f %3.1f %3.1f\n",trackloc[0],trackloc[1],trackloc[2]);
273
274
275     iChamber->LocaltoGlobal(trackloc,trackglob);
276        
277     //printf("Transformation 2: %3.1f %3.1f %3.1f\n",trackglob[0],trackglob[1],trackglob[2]);
278     
279     cx=trackloc[0];
280     cy=trackloc[2];
281      
282
283     TClonesArray *pDigits = pRICH->DigitsAddress(nch-1);   
284     Int_t ndigits = pDigits->GetEntriesFast();
285     
286     //printf("Got %d digits\n",ndigits);
287
288     counter=0;
289     printf("Starting calculations\n");
290     for(Float_t theta=0;theta<kMaxTheta;theta+=steptheta)
291       {         
292         //printf(".");
293         for(Float_t phi=0;phi<=kMaxPhi;phi+=stepphi)
294           {             
295             //printf("Phi:%3.1f\n", phi*180/kPi);
296             counter1=0;
297             for (Int_t dig=0;dig<ndigits;dig++)
298               { 
299                 points=(AliRICHDigit*) pDigits->UncheckedAt(dig);
300                 segmentation->GetPadC(points->fPadX, points->fPadY,x, y, z);
301                 x=x-cx;
302                 y=y-cy;
303                 radius=TMath::Sqrt(TMath::Power(x,2)+TMath::Power(y,2));
304
305                 if(radius>4)
306                   {
307                     //if(theta==0 && phi==0)
308                       //{
309                         //printf("Radius: %f, Max Radius: %f\n",radius,kCorr*kHeight*tan(theta+kMaxOmega)*3/4);
310                         meanradius+=radius;
311                         counter++;
312                       //}
313                     
314                     if (radius<2*kHeight*tan(theta+kMaxOmega)*3/4)
315                       {
316                         
317                         if(phi==0)
318                           {
319                             //printf("Radius: %f, Max Radius: %f\n",radius,2*kHeight*tan(theta+kMaxOmega)*3/4);
320                             //printf("Loaded digit %d with coordinates x:%f, y%f\n",dig,x,y);
321                             //printf("Using digit %d, for theta %f\n",dig,theta);
322                           }
323                         
324                         counter1++;
325
326                         l=kHeight/cos(theta);
327                         
328                         //x=x*kCorr;
329                         //y=y*kCorr;
330                         /*if(SnellAngle(theta+omega)<999)
331                           {
332                             //printf("(Angle)/(Snell angle):%f\n",(theta+omega)/SnellAngle(theta+omega));
333                             x=x*(theta+omega)/SnellAngle(theta+omega);
334                             y=y*(theta+omega)/SnellAngle(theta+omega);
335                           }
336                         else
337                           {
338                             x=0;
339                             y=0;
340                           }*/
341
342                         //main calculation
343
344                         DigitsXY->Fill(x,y,(float) 1);
345
346                         theta1=SnellAngle(theta)*1.5;
347                 
348                         aux1=-y*sin(phi)+x*cos(phi);
349                         aux2=y*cos(phi)+x*sin(phi);
350                         aux3=( TMath::Power(aux1,2)+TMath::Power(cos(theta1)*aux2 ,2))/TMath::Power(sin(theta1)*aux2+l,2);
351                         omega=atan(sqrt(aux3));
352                         
353                         //omega is distorted, theta1 is distorted
354
355                         if(InvSnellAngle(theta+omega)<999)
356                           {
357                             omega1=InvSnellAngle(omega+theta1) - theta;
358                             //theta1=InvSnellAngle(omega+theta) - omega1;
359                             //omega1=kCorr*omega;
360                             
361                             kCorr=InvSnellAngle(omega+theta)/(omega+theta);
362                             theta1=kCorr*theta/1.4;
363                             //if(phi==0)
364                               //printf("Omega:%f Theta:%f Omega1:%f Theta1:%f ISA(o+t):%f ISA(t):%f\n",omega*180/kPi,theta*180/kPi,omega1*180/kPi,theta1*180/kPi,InvSnellAngle(omega+theta)*180/kPi,InvSnellAngle(theta)*180/kPi);
365                           }
366                         else
367                           {
368                             omega1=0;
369                             theta1=0;
370                           }
371                         
372                         //printf("Omega:%f\n",omega);
373
374
375                         //if(SnellAngle(theta+omega)<999)
376                           //printf("(Angle)/(Snell angle):%f\n",(theta+omega)/SnellAngle(theta+omega));
377                         if(theta==0 && phi==0)
378                           {
379                             //printf("Omega: %f Corrected Omega: %f\n",omega, omega/kCorr);
380                             //omega=omega/kCorr;
381                           }
382                         
383                         //cout<<"\ni="<<i<<" theta="<<Int_t(2*theta*dimension/kPi)<<" phi="<<Int_t(2*phi*dimension/kPi)<<" omega="<<Int_t(2*omega*dimension/kPi)<<endl<<endl;
384                         //{Int_t lixo;cin>>lixo;}
385                         if(omega1<kMaxOmega && omega1>kMinOmega)
386                           {
387                             //printf("Omega found:%f\n",omega);
388                             omega1=omega1-kMinOmega;
389                             
390                             //printf("Omega: %f Theta: %3.1f Phi:%3.1f\n",omega, theta*180/kPi, phi*180/kPi);
391
392                             bintheta=theta1*kDimensionTheta/kMaxTheta;
393                             binphi=phi*kDimensionPhi/kMaxPhi;
394                             binomega=omega1*kDimensionOmega/(kMaxOmega-kMinOmega);
395
396                             if(Int_t(bintheta+0.5)==Int_t(bintheta))
397                               inttheta=Int_t(bintheta);
398                             else
399                               inttheta=Int_t(bintheta+0.5);
400
401                             if(Int_t(binomega+0.5)==Int_t(binomega))
402                               intomega=Int_t(binomega);
403                             else
404                               intomega=Int_t(binomega+0.5);
405                             
406                             if(Int_t(binphi+0.5)==Int_t(binphi))
407                               intphi=Int_t(binphi);
408                             else
409                               intphi=Int_t(binphi+0.5);
410                                                  
411                             //printf("Point added at %d %d %d\n",inttheta,intphi,intomega);
412                             point[inttheta][intphi][intomega]+=1;
413                             //printf("Omega stored:%d\n",intomega);
414                             Points->Fill(inttheta,intphi,intomega,(float) 1);
415                             ThetaPhi->Fill(inttheta,intphi,(float) 1);
416                             OmegaTheta->Fill(inttheta,intomega,(float) 1);
417                             OmegaPhi->Fill(intphi,intomega,(float) 1);
418                             //printf("Filling at %d %d %d\n",Int_t(theta*kDimensionTheta/kMaxTheta),Int_t(phi*kDimensionPhi/kMaxPhi),Int_t(omega*kDimensionOmega/kMaxOmega));
419                           }
420                         //if(omega<kMaxOmega)point[Int_t(theta)][Int_t(phi)][Int_t(omega)]+=1;
421                       }
422                   }
423               }
424           }
425         //printf("Used %d digits for theta %3.1f\n",counter1, theta*180/kPi);
426       }
427
428     meanradius=meanradius/counter;
429     printf("Mean radius:%f, counter:%d\n",meanradius,counter);
430     rechit[5]=meanradius;
431     printf("Used %d digits\n",counter1);
432     //printf("\n");
433
434     if(nev<20)
435       {
436         if(nev==0)
437           {
438             fc1->cd(1);
439             Points->Draw();
440             fc1->cd(2);
441             ThetaPhi->Draw();
442             fc1->cd(3);
443             OmegaTheta->Draw();
444             fc1->cd(4);
445             OmegaPhi->Draw();
446             fc3->cd();
447             DigitsXY->Draw();
448           }
449         else
450           {
451             //fc1->cd(1);
452             //Points->Draw("same");
453             //fc1->cd(2);
454             //ThetaPhi->Draw("same");
455             //fc1->cd(3);
456             //OmegaTheta->Draw("same");
457             //fc1->cd(4);
458             //OmegaPhi->Draw("same");   
459           }
460       }
461         
462     
463     //SPOT execute twice
464     for(Int_t s=0;s<kSpot;s++)
465       {
466         printf("     Applying Spot algorithm, pass %d\n", s);
467         
468         //buffer copy
469         for(i=0;i<=kDimensionTheta;i++)
470           {
471             for(j=0;j<=kDimensionPhi;j++)
472               {
473                 for(k=0;k<=kDimensionOmega;k++)
474                   {
475                     point1[i][j][k]=point[i][j][k];     
476                   }
477               }
478           }
479
480         //SPOT algorithm                        
481         for(i=1;i<kDimensionTheta-1;i++)
482           {
483             for(j=1;j<kDimensionPhi-1;j++)
484               {
485                 for(k=1;k<kDimensionOmega-1;k++)
486                   {
487                     if((point[i][k][j]>point[i-1][k][j])&&(point[i][k][j]>point[i+1][k][j])&&
488                        (point[i][k][j]>point[i][k-1][j])&&(point[i][k][j]>point[i][k+1][j])&&
489                        (point[i][k][j]>point[i][k][j-1])&&(point[i][k][j]>point[i][k][j+1]))
490                       {
491                         //cout<<"SPOT"<<endl;
492                         //Execute SPOT on point                                                                                         
493                         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]));    
494                         point1[i-1][k][j]=Int_t(SPOTp*point[i-1][k][j]);
495                         point1[i+1][k][j]=Int_t(SPOTp*point[i+1][k][j]);
496                         point1[i][k-1][j]=Int_t(SPOTp*point[i][k-1][j]);
497                         point1[i][k+1][j]=Int_t(SPOTp*point[i][k+1][j]);
498                         point1[i][k][j-1]=Int_t(SPOTp*point[i][k][j-1]);
499                         point1[i][k][j+1]=Int_t(SPOTp*point[i][k][j+1]);
500                       }
501                   }
502               }
503           }
504         
505         //copy from buffer copy
506         counter1=0;
507         for(i=1;i<kDimensionTheta;i++)
508           {
509             for(j=1;j<kDimensionPhi;j++)
510               {
511                 for(k=1;k<kDimensionOmega;k++)
512                   {
513                     point[i][j][k]=point1[i][j][k];
514                     if(nev<20)
515                       {
516                         if(s==kSpot-1)
517                           {
518                             if(point1[i][j][k] != 0)
519                               {
520                                 SpotPoints->Fill(i,j,k,(float) point1[i][j][k]);
521                                 //printf("Random number %f\n",random->Rndm(2));
522                                 //if(random->Rndm() < .2)
523                                   //{
524                                 SpotThetaPhi->Fill(i,j,(float) point1[i][j][k]);
525                                 SpotOmegaTheta->Fill(i,k,(float) point1[i][j][k]);
526                                 SpotOmegaPhi->Fill(j,k,(float) point1[i][j][k]);
527                                     counter1++;
528                                   //}
529                                 //printf("Filling at %d %d %d value %f\n",i,j,k,(float) point1[i][j][k]);
530                               }
531                           }
532                       }
533                     //if(point1[i][j][k] != 0)
534                       //printf("Last transfer point: %d, point1, %d\n",point[i][j][k],point1[i][j][k]);
535                   }
536               }
537           }
538       }
539     
540     //printf("Filled %d cells\n",counter1);
541
542     if(nev<20)
543       {
544         if(nev==0)
545           {
546             fc2->cd(1);
547             SpotPoints->Draw();
548             fc2->cd(2);
549             SpotThetaPhi->Draw();
550             fc2->cd(3);
551             SpotOmegaTheta->Draw();
552             fc2->cd(4);
553             SpotOmegaPhi->Draw();
554           }
555         else
556           {
557             //fc2->cd(1);
558             //SpotPoints->Draw("same");
559             //fc2->cd(2);
560             //SpotThetaPhi->Draw("same");
561             //fc2->cd(3);
562             //SpotOmegaTheta->Draw("same");
563             //fc2->cd(4);
564             //SpotOmegaPhi->Draw("same");       
565           }
566       }
567     
568     
569     //Identification is equivalent to maximum determination
570     max=0;maxi=0;maxj=0;maxk=0;
571     
572     printf("     Proceeding to identification");
573     
574     for(i=0;i<kDimensionTheta;i++)
575       for(j=0;j<kDimensionPhi;j++)
576         for(k=0;k<kDimensionOmega;k++)
577             if(point[i][j][k]>max)
578               {
579                 //cout<<"maxi="<<i*90/dimension<<" maxj="<<j*90/dimension<<" maxk="<<k*kMaxOmega/dimension*180/kPi<<" max="<<max<<endl;
580                 maxi=i;maxj=j;maxk=k;
581                 max=point[i][j][k];
582                 printf(".");
583                 //printf("Max Omega %d, Max Theta %d, Max Phi %d (%d counts)\n",maxk,maxi,maxj,max);
584               }
585     printf("\n");
586     
587     Float_t FinalOmega = maxk*(kMaxOmega-kMinOmega)/kDimensionOmega; 
588     Float_t FinalTheta = maxi*kMaxTheta/kDimensionTheta;
589     Float_t FinalPhi = maxj*kMaxPhi/kDimensionPhi;
590
591     FinalOmega += kMinOmega;
592     
593     //printf("Detected angle for height %3.1f and for center %3.1f %3.1f:%f\n",h,cx,cy,maxk*kPi/(kDimensionTheta*4));
594     printf("     Indentified angles: cerenkov - %f, theta - %3.1f, phi - %3.1f (%f activation)\n", FinalOmega, FinalTheta*180/kPi, FinalPhi*180/kPi, max);
595     //printf("Detected angle for height %3.1f and for center %3.1f %3.1f:%f\n",kHeight,cx,cy,maxk);
596
597     //fscanf(omegas,"%f",&realomega);
598     //fscanf(thetas,"%f",&realtheta);
599     //printf("Real Omega: %f",realomega);                       
600     //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;           
601     
602     //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));
603     
604     /*for(j=0;j<np;j++)
605       pointpp(maxj*90/kDimensionTheta,maxi*90/kDimensionPhi,maxk*kMaxOmega/kDimensionOmega*180/kPi,cx,cy);//Generates a point on the elipse*/               
606
607
608     //Start filling rec. hits
609     
610     rechit[0] = FinalTheta;
611     rechit[1] = 90*kPi/180 + FinalPhi;
612     rechit[2] = FinalOmega;
613     rechit[3] = cx;
614     rechit[4] = cy;
615
616     //CreatePoints(FinalTheta, 270*kPi/180 + FinalPhi, FinalOmega, kHeight);
617        
618     //printf ("track %d, theta %f, phi %f, omega %f\n\n\n",track,rechit[0],rechit[1],rechit[2]);
619     
620     // fill rechits
621     pRICH->AddRecHit3D(nch-1,rechit);
622     //printf("rechit %f %f %f %f %f\n",rechit[0],rechit[1],rechit[2],rechit[3],rechit[4]);
623     //printf("Chamber:%d",nch);
624   }                     
625   //printf("\n\n\n\n");
626   gAlice->TreeR()->Fill();
627   TClonesArray *fRec;
628   for (i=0;i<kNCH;i++) {
629     fRec=pRICH->RecHitsAddress3D(i);
630     int ndig=fRec->GetEntriesFast();
631     printf ("Chamber %d, rings %d\n",i+1,ndig);
632   }
633   pRICH->ResetRecHits3D();
634
635   free_i3tensor(point,0,kDimensionTheta,0,kDimensionPhi,0,kDimensionOmega);
636   free_i3tensor(point1,0,kDimensionTheta,0,kDimensionPhi,0,kDimensionOmega);
637 }
638
639
640
641 Float_t AliRICHDetect:: Area(Float_t theta,Float_t omega)
642 {
643
644 //
645 // Calculates area of an ellipse for given incidence angles    
646
647
648     Float_t area;
649     const Float_t kHeight=9.25;                       //Distance from Radiator to Pads in pads
650     
651     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);
652     
653     return (area);
654 }
655
656
657 Int_t ***AliRICHDetect::i3tensor(long nrl, long nrh, long ncl, long nch, long ndl, long ndh)
658 // allocate a Int_t 3tensor with range t[nrl..nrh][ncl..nch][ndl..ndh] 
659 {
660   long i,j,nrow=nrh-nrl+1,ncol=nch-ncl+1,ndep=ndh-ndl+1;
661   Int_t ***t;
662   
663   int NR_END=1; 
664
665   // allocate pointers to pointers to rows 
666   t=(Int_t ***) malloc((size_t)((nrow+NR_END)*sizeof(Int_t**)));
667   if (!t) printf("allocation failure 1 in f3tensor()");
668   t += NR_END;
669   t -= nrl;
670   
671   // allocate pointers to rows and set pointers to them 
672   t[nrl]=(Int_t **) malloc((size_t)((nrow*ncol+NR_END)*sizeof(Int_t*)));
673   if (!t[nrl]) printf("allocation failure 2 in f3tensor()");
674   t[nrl] += NR_END;
675   t[nrl] -= ncl;
676   
677   // allocate rows and set pointers to them 
678   t[nrl][ncl]=(Int_t *) malloc((size_t)((nrow*ncol*ndep+NR_END)*sizeof(Int_t)));
679   if (!t[nrl][ncl]) printf("allocation failure 3 in f3tensor()");
680   t[nrl][ncl] += NR_END;
681   t[nrl][ncl] -= ndl;
682   
683   for(j=ncl+1;j<=nch;j++) t[nrl][j]=t[nrl][j-1]+ndep;
684   for(i=nrl+1;i<=nrh;i++) {
685     t[i]=t[i-1]+ncol;
686     t[i][ncl]=t[i-1][ncl]+ncol*ndep;
687     for(j=ncl+1;j<=nch;j++) t[i][j]=t[i][j-1]+ndep;
688   }
689   
690   // return pointer to array of pointers to rows 
691   return t;
692 }
693
694 void AliRICHDetect::free_i3tensor(int ***t, long nrl, long nrh, long ncl, long nch,long ndl, long ndh)
695 // free a Int_t f3tensor allocated by i3tensor()
696 {
697   int NR_END=1; 
698
699   free((char*) (t[nrl][ncl]+ndl-NR_END));
700   free((char*) (t[nrl]+ncl-NR_END));
701   free((char*) (t+nrl-NR_END));
702 }
703
704
705 Float_t AliRICHDetect:: SnellAngle(Float_t iangle)
706
707
708 // Compute the Snell angle
709
710   Float_t nfreon  = 1.295;
711   Float_t nquartz = 1.585;
712   Float_t ngas = 1;
713
714   Float_t sinrangle;
715   Float_t rangle;
716   Float_t a1, a2;
717
718   a1=nfreon/nquartz;
719   a2=nquartz/ngas;
720
721   sinrangle = a1*a2*sin(iangle);
722
723   if(sinrangle>1.) {
724      rangle = 999.;
725      return rangle;
726   }
727   
728   rangle = asin(sinrangle);  
729   return rangle;
730 }
731
732 Float_t AliRICHDetect:: InvSnellAngle(Float_t rangle)
733
734
735 // Compute the inverse Snell angle
736
737   Float_t nfreon  = 1.295;
738   Float_t nquartz = 1.585;
739   Float_t ngas = 1;
740
741   Float_t siniangle;
742   Float_t iangle;
743   Float_t a1,a2;
744
745   a1=nfreon/nquartz;
746   a2=nquartz/ngas;
747
748   siniangle = sin(rangle)/(a1*a2);
749   iangle = asin(siniangle);
750
751   if(siniangle>1.) {
752      iangle = 999.;
753      return iangle;
754   }
755   
756   iangle = asin(siniangle);
757   return iangle;
758 }
759
760
761
762 //________________________________________________________________________________
763 void AliRICHDetect::CreatePoints(Float_t theta, Float_t phi, Float_t omega, Float_t h)
764 {
765   
766   // Create points along the ellipse equation
767   
768   Int_t s1,s2;
769   Float_t fiducial=h*TMath::Tan(omega+theta), l=h/TMath::Cos(theta), xtrial, y=0, c0, c1, c2;
770   //TRandom *random=new TRandom();
771   
772   static TH2F *REllipse = new TH2F("REllipse","Reconstructed ellipses",150,-25,25,150,-25,25);
773
774   for(Float_t i=0;i<1000;i++)
775     {
776       
777       Float_t counter=0;
778       
779       c0=0;c1=0;c2=0;
780       while((c1*c1-4*c2*c0)<=0 && counter<1000)
781         {
782           //Choose which side to go...
783           if(i>250 && i<750) s1=1; 
784           //if (gRandom->Rndm(1)>.5) s1=1;
785           else s1=-1;
786           //printf("s1:%d\n",s1);
787           //Trial a y
788           y=s1*i*gRandom->Rndm(Int_t(fiducial/50));
789           //printf("Fiducial %f  for omega:%f theta:%f phi:%f\n",fiducial,omega,theta,fphi);
790           Float_t alfa1=theta;
791           Float_t theta1=phi;
792           Float_t omega1=omega;
793           
794           //Solve the eq for a trial x
795           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);
796           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);
797           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);
798           //cout<<"Trial: y="<<y<<"c0="<<c0<<" c1="<<c1<<" c2="<<c2<<endl;
799           //printf("Result:%f\n\n",c1*c1-4*c2*c0);
800           //i+=.01;
801           counter +=1;
802         }
803       
804       if (counter>=1000)
805         y=0; 
806
807       //Choose which side to go...
808       //if(gRandom->Rndm(1)>.5) s=1; 
809       //else s=-1;
810       if(i>500) s2=1;
811       //if (gRandom->Rndm(1)>.5) s2=1;
812       else s2=-1;
813       xtrial=(-c1+s2*TMath::Sqrt(c1*c1-4*c2*c0))/(2*c2);
814       //cout<<"x="<<xtrial<<" y="<<cy+y<<endl;
815       //printf("Coordinates: %f %f\n",xtrial,fCy+y);
816
817       REllipse->Fill(xtrial,y);
818
819       //printf("Coordinates: %f %f %f\n",vectorGlob[0],vectorGlob[1],vectorGlob[2]);
820     }
821   
822   fc3->cd(2);
823   REllipse->Draw();
824 }