]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RICH/AliRICHDetect.cxx
Container for calibration data
[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.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)
20
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
23
24   Revision 1.14  2001/05/14 13:25:54  hristov
25   stdlib.h included (for Alpha)
26
27   Revision 1.13  2001/05/10 12:26:31  jbarbosa
28   Totally reworked version of reconstruction algorithm.
29
30   Revision 1.12  2001/02/27 22:15:03  jbarbosa
31   Removed compiler warning.
32
33   Revision 1.11  2001/02/27 15:21:46  jbarbosa
34   Transition to SDigits.
35
36   Revision 1.10  2001/02/13 20:39:06  jbarbosa
37   Changes to make it work with new IO.
38
39   Revision 1.9  2001/01/22 21:39:11  jbarbosa
40   Several tune-ups
41
42   Revision 1.8  2000/11/15 15:52:53  jbarbosa
43   Turned on spot algorithm.
44
45   Revision 1.7  2000/11/01 15:37:05  jbarbosa
46   Updated to use its own rec. point object.
47
48   Revision 1.6  2000/10/02 21:28:12  fca
49   Removal of useless dependecies via forward declarations
50
51   Revision 1.5  2000/06/30 16:30:28  dibari
52   Disabled writing to rechits.
53
54   Revision 1.4  2000/06/15 15:46:59  jbarbosa
55   Corrected compilation errors on HP-UX (replaced pow with TMath::Power)
56
57   Revision 1.3  2000/06/13 13:15:41  jbarbosa
58   Still some code cleanup done (variable names)
59
60   Revision 1.2  2000/06/12 15:19:30  jbarbosa
61   Cleaned up version.
62
63   Revision 1.1  2000/04/19 13:05:14  morsch
64   J. Barbosa's spot reconstruction algorithm.
65
66 */
67
68 #include <stdlib.h>
69
70
71 #include "AliRICH.h"
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"
79 #include "AliRun.h"
80 #include "TParticle.h"
81 #include "TTree.h"
82 #include "TMath.h"
83 #include "TRandom.h"
84 #include "TH3.h"
85 #include "TH2.h"
86 #include "TCanvas.h"
87 #include <TStyle.h>
88
89
90 #include "malloc.h"
91
92
93 ClassImp(AliRICHDetect)
94 //___________________________________________
95 AliRICHDetect::AliRICHDetect() : TObject()
96 {
97
98 // Default constructor 
99
100   fc1 = 0;
101   fc2 = 0;
102   fc3 = 0;
103
104 }
105
106 //___________________________________________
107 AliRICHDetect::AliRICHDetect(const char *name, const char *title)
108     : TObject()
109 {
110
111   TStyle *mystyle=new TStyle("Plain","mystyle");
112   mystyle->SetPalette(1,0);
113   mystyle->cd();
114
115
116   fc1= new TCanvas("c1","Reconstructed points",50,50,300,350);
117   fc1->Divide(2,2);
118   fc2= new TCanvas("c2","Reconstructed points after SPOT",370,50,300,350);
119   fc2->Divide(2,2); 
120   fc3= new TCanvas("c3","Used Digits",690,50,300,350);
121   fc4= new TCanvas("c4","Mesh activation data",50,430,600,350);
122   fc4->Divide(2,1);
123  
124
125 }
126
127 //___________________________________________
128 AliRICHDetect::~AliRICHDetect()
129 {
130     
131 // Destructor
132
133 }
134
135
136 void AliRICHDetect::Detect(Int_t nev, Int_t type)
137 {       
138     
139 //
140 // Detection algorithm
141
142
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;
149   Int_t i,j,k;
150
151   AliRICH *pRICH  = (AliRICH*)gAlice->GetDetector("RICH");
152   AliRICHSegmentationV0*  segmentation;
153   AliRICHChamber*       iChamber;
154   AliRICHGeometry*  geometry;
155   
156   iChamber = &(pRICH->Chamber(0));
157   segmentation=(AliRICHSegmentationV0*) iChamber->GetSegmentationModel(0);
158   geometry=iChamber->GetGeometryModel();
159  
160   
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)
163   
164   const Float_t kPi=TMath::Pi();                
165   
166   const Float_t kHeight=geometry->GetRadiatorToPads(); //Distance from Radiator to Pads in centimeters
167   //printf("Distance to Pads:%f\n",kHeight);
168  
169   const Int_t kSpot=0;                                 //number of passes with spot algorithm
170   
171   const Int_t kDimensionTheta=100;                     //Matrix dimension for angle Detection
172   const Int_t kDimensionPhi=100;
173   const Int_t kDimensionOmega=100;
174   
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;
182
183   Float_t rechit[6];                                 //Reconstructed point data
184
185   Int_t ***point = i3tensor(0,kDimensionTheta,0,kDimensionPhi,0,kDimensionOmega);
186   Int_t ***point1 = i3tensor(0,kDimensionTheta,0,kDimensionPhi,0,kDimensionOmega);
187   
188   steptheta=(kMaxTheta-kMinTheta)/kDimensionTheta;
189   stepphi=(kMaxPhi-kMinPhi)/kDimensionPhi;
190
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");
225
226   Int_t ntracks = (Int_t)gAlice->TreeH()->GetEntries();
227    
228   Float_t trackglob[3];
229   Float_t trackloc[3];
230
231   //printf("Area de uma elipse com teta 0 e Omega 45:%f",Area(0,45));
232    
233   Int_t track;
234         
235   for (track=0; track<ntracks;track++) {
236     gAlice->ResetHits();
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;
244     //Int_t npoints=0;
245     
246     Int_t counter=0, counter1=0;
247     //Initialization
248     for(i=0;i<kDimensionTheta;i++)
249       {
250         for(j=0;j<kDimensionPhi;j++)
251           {
252             for(k=0;k<kDimensionOmega;k++)
253               {
254                 counter++;
255                 point[i][j][k]=0;
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)));
260               }
261           }
262       }
263
264     Int_t ncerenkovs = pRICH->Cerenkovs()->GetEntriesFast();
265     
266     
267     originalOmega = 0;
268     counter = 0;
269     
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)
276         {
277           counter +=1;
278           originalOmega += cherenkov;
279           //printf("%f\n",cherenkov);
280         }
281     }
282     
283     printf("Cerenkovs       : %d\n",counter);
284     
285     if(counter != 0)    //if there are cerenkovs
286       {
287         originalOmega = originalOmega/counter;
288         printf("Original omega: %f\n",originalOmega);
289      
290
291
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();
299         
300         printf("\n--------------------------------------\n");
301         printf("Chamber %d, track %d\n", nch, track);
302
303         
304         iChamber = &(pRICH->Chamber(nch-1));
305     
306         //printf("Nch:%d\n",nch);
307
308         iChamber->GlobaltoLocal(trackglob,trackloc);
309     
310         iChamber->LocaltoGlobal(trackloc,trackglob);
311        
312        
313         cx=trackloc[0];
314         cy=trackloc[2];
315      
316         AliRICHDigit *points = 0;
317         TClonesArray *pDigits = pRICH->DigitsAddress(nch-1);   
318         
319         AliRICHRawCluster *cluster =0;
320         TClonesArray *pClusters = pRICH->RawClustAddress(nch-1);
321
322         Int_t maxcycle=0;
323
324         //digitize from digits
325         
326         if(type==0)
327           {
328             gAlice->TreeD()->GetEvent(0);
329             Int_t ndigits = pDigits->GetEntriesFast();
330             maxcycle=ndigits;
331             printf("Got %d digits\n",ndigits);
332           }
333         
334         //digitize from clusters
335         
336         if(type==1)
337           {
338             Int_t nent=(Int_t)gAlice->TreeR()->GetEntries();
339             gAlice->TreeR()->GetEvent(nent-1);
340             Int_t nclusters = pClusters->GetEntriesFast();
341             maxcycle=nclusters;
342             printf("Got %d clusters\n",nclusters);
343           }
344
345
346
347         
348         counter=0;
349         printf("Starting calculations\n");
350         printf("           Start                                                                                              Finish\n");
351         printf("Progress:     ");
352         for(Float_t theta=0;theta<kMaxTheta;theta+=steptheta)
353           {             
354             printf(".");
355             for(Float_t phi=0;phi<=kMaxPhi;phi+=stepphi)
356               {         
357                 //printf("Phi:%3.1f\n", phi*180/kPi);
358                 counter1=0;
359                 for (Int_t cycle=0;cycle<maxcycle;cycle++)
360                   {     
361                     
362                     if(type==0)
363                       {
364                         points=(AliRICHDigit*) pDigits->UncheckedAt(cycle);
365                         segmentation->GetPadC(points->PadX(), points->PadY(),x, y, z);
366                       }
367                     
368                     if(type==1)
369                       {
370                         cluster=(AliRICHRawCluster*) pClusters->UncheckedAt(cycle);
371                         x=cluster->fX;
372                         y=cluster->fY;
373                         q=cluster->fQ;
374                       }
375                     
376                     if(type ==0 || q > 100)
377                       
378                       {
379                         
380                         x=x-cx;
381                         y=y-cy;
382                         radius=TMath::Sqrt(TMath::Power(x,2)+TMath::Power(y,2));
383                         
384                         //calculation of relative phi angle of digit
385                         
386                         phi_relative = acos(y/radius);
387                         phi_relative = TMath::Abs(phi_relative - phi);
388                         
389                         
390                         if(radius>4)
391                           {
392                             meanradius+=radius;
393                             counter++;
394                             
395                             
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))
399                               {
400                                 
401                                 if(phi==0)
402                                   {
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);
406                                   }
407                                 
408                                 counter1++;
409                                 
410                                 l=kHeight/cos(theta);
411                                 
412                                 //main calculation
413                                 
414                                 DigitsXY->Fill(x,y,(float) 1);
415                                 
416                                 theta1=SnellAngle(theta);
417                                 
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));
422                                 
423                                 //omega is distorted, theta1 is distorted
424                                 
425                                 if(InvSnellAngle(omega+TMath::Abs(cos(phi_relative))*theta1)<999)
426                                   {
427                                     omega1=InvSnellAngle(omega+TMath::Abs(cos(phi_relative))*theta);
428                                     theta1=InvSnellAngle(theta1);
429                                     
430                                   }
431                                 else
432                                   {
433                                     omega1=0;
434                                     theta1=0;
435                                   }
436                                 
437                                 if(omega1<kMaxOmega && omega1>kMinOmega)
438                                   {
439                                 //printf("Omega found:%f\n",omega);
440                                     omega1=omega1-kMinOmega;
441                                     
442                                 //printf("Omega: %f Theta: %3.1f Phi:%3.1f\n",omega, theta*180/kPi, phi*180/kPi);
443                                     
444                                     bintheta=theta1*kDimensionTheta/kMaxTheta;
445                                     binphi=phi*kDimensionPhi/kMaxPhi;
446                                     binomega=omega1*kDimensionOmega/(kMaxOmega-kMinOmega);
447                                     
448                                     if(Int_t(bintheta+0.5)==Int_t(bintheta))
449                                       inttheta=Int_t(bintheta);
450                                     else
451                                       inttheta=Int_t(bintheta+0.5);
452                                     
453                                     if(Int_t(binomega+0.5)==Int_t(binomega))
454                                       intomega=Int_t(binomega);
455                                     else
456                                       intomega=Int_t(binomega+0.5);
457                                     
458                                     if(Int_t(binphi+0.5)==Int_t(binphi))
459                                       intphi=Int_t(binphi);
460                                     else
461                                       intphi=Int_t(binphi+0.5);
462                                     
463                                 //printf("Point added at %d %d %d\n",inttheta,intphi,intomega);
464                                     
465                                     if(type==0)
466                                       point[inttheta][intphi][intomega]+=1;
467                                     if(type==1)
468                                       point[inttheta][intphi][intomega]+=(Int_t)(q);
469                                     
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));
476                                   }
477                                 //if(omega<kMaxOmega)point[Int_t(theta)][Int_t(phi)][Int_t(omega)]+=1;
478                               }
479                           }
480                       }
481                     
482                   }
483                 //printf("Used %d digits for theta %3.1f\n",counter1, theta*180/kPi);
484               }
485             
486           }
487
488         printf("\n");
489
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);
494         //printf("\n");
495
496         if(nev<2)
497           {
498             if(nev==0)
499               {
500                 fc1->cd(1);
501                 Points->Draw("colz");
502                 fc1->cd(2);
503                 ThetaPhi->Draw("colz");
504                 fc1->cd(3);
505                 OmegaTheta->Draw("colz");
506                 fc1->cd(4);
507                 OmegaPhi->Draw("colz");
508                 fc3->cd();
509                 DigitsXY->Draw("colz");
510               }
511             else
512               {
513                 //fc1->cd(1);
514                 //Points->Draw("same");
515                 //fc1->cd(2);
516                 //ThetaPhi->Draw("same");
517                 //fc1->cd(3);
518                 //OmegaTheta->Draw("same");
519                 //fc1->cd(4);
520                 //OmegaPhi->Draw("same");       
521               }
522           }
523         
524     
525         //SPOT execute twice
526         for(Int_t s=0;s<kSpot;s++)
527           {
528             printf("     Applying Spot algorithm, pass %d\n", s);
529         
530         //buffer copy
531         for(i=0;i<=kDimensionTheta;i++)
532           {
533             for(j=0;j<=kDimensionPhi;j++)
534               {
535                 for(k=0;k<=kDimensionOmega;k++)
536                   {
537                     point1[i][j][k]=point[i][j][k];     
538                   }
539               }
540           }
541
542         //SPOT algorithm                        
543         for(i=1;i<kDimensionTheta-1;i++)
544           {
545             for(j=1;j<kDimensionPhi-1;j++)
546               {
547                 for(k=1;k<kDimensionOmega-1;k++)
548                   {
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]))
552                       {
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]);
562                       }
563                   }
564               }
565           }
566         
567         //copy from buffer copy
568         counter1=0;
569         for(i=1;i<kDimensionTheta;i++)
570           {
571             for(j=1;j<kDimensionPhi;j++)
572               {
573                 for(k=1;k<kDimensionOmega;k++)
574                   {
575                     point[i][j][k]=point1[i][j][k];
576                     if(nev<20)
577                       {
578                         if(s==kSpot-1)
579                           {
580                             if(point1[i][j][k] != 0)
581                               {
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)
585                                   //{
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]);
589                                     counter1++;
590                                   //}
591                                 //printf("Filling at %d %d %d value %f\n",i,j,k,(float) point1[i][j][k]);
592                               }
593                           }
594                       }
595                     //if(point1[i][j][k] != 0)
596                       //printf("Last transfer point: %d, point1, %d\n",point[i][j][k],point1[i][j][k]);
597                   }
598               }
599           }
600       }
601     
602     //printf("Filled %d cells\n",counter1);
603
604         if(nev<2)
605           {
606             if(nev==0)
607               {
608                 fc2->cd(1);
609                 SpotPoints->Draw("colz");
610                 fc2->cd(2);
611                 SpotThetaPhi->Draw("colz");
612                 fc2->cd(3);
613                 SpotOmegaTheta->Draw("colz");
614                 fc2->cd(4);
615                 SpotOmegaPhi->Draw("colz");
616               }
617             else
618               {
619                 //fc2->cd(1);
620                 //SpotPoints->Draw("same");
621                 //fc2->cd(2);
622                 //SpotThetaPhi->Draw("same");
623                 //fc2->cd(3);
624                 //SpotOmegaTheta->Draw("same");
625                 //fc2->cd(4);
626                 //SpotOmegaPhi->Draw("same");   
627               }
628           }
629     
630     
631         //Identification is equivalent to maximum determination
632         max=0;maxi=0;maxj=0;maxk=0;
633         
634         printf("     Proceeding to identification");
635         
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)
640                 {
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;
643                 max=point[i][j][k];
644                 printf(".");
645                 //printf("Max Omega %d, Max Theta %d, Max Phi %d (%d counts)\n",maxk,maxi,maxj,max);
646                 }
647         printf("\n");
648     
649         Float_t FinalOmega = maxk*(kMaxOmega-kMinOmega)/kDimensionOmega; 
650         Float_t FinalTheta = maxi*kMaxTheta/kDimensionTheta;
651         Float_t FinalPhi = maxj*kMaxPhi/kDimensionPhi;
652
653         FinalOmega += kMinOmega;
654     
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);
658
659         AngleAct->Fill(FinalOmega, (float) max);
660         Activation->Fill(max, (float) 1);
661
662         if(nev==0)
663               {
664                 fc4->cd(1);
665                 AngleAct->Draw();
666                 fc4->cd(2);
667                 Activation->Draw();
668               }
669             else
670               {
671                 fc4->cd(1);
672                 AngleAct->Draw("same");
673                 fc4->cd(2);
674                 Activation->Draw("same");
675               }
676           
677
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;               
682     
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));
684     
685     /*for(j=0;j<np;j++)
686       pointpp(maxj*90/kDimensionTheta,maxi*90/kDimensionPhi,maxk*kMaxOmega/kDimensionOmega*180/kPi,cx,cy);//Generates a point on the elipse*/               
687
688
689     //Start filling rec. hits
690     
691     rechit[0] = FinalTheta;
692     rechit[1] = FinalPhi - 90*kPi/180;
693     rechit[2] = FinalOmega;
694     rechit[3] = cx;
695     rechit[4] = cy;
696
697     //CreatePoints(FinalTheta, 270*kPi/180 + FinalPhi, FinalOmega, kHeight);
698        
699     printf ("track %d, theta %f, phi %f, omega %f\n\n\n",track,rechit[0],rechit[1]*180/kPi,rechit[2]);
700     
701     // fill rechits
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);
705       }
706
707     else   //if no cerenkovs
708       
709       {
710         
711         rechit[0] = 0;
712         rechit[1] = 0;
713         rechit[2] = 0;
714         rechit[3] = 0;
715         rechit[4] = 0;
716
717       }
718
719   }
720
721   if(type==1)  //reco from clusters
722     {
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);
728     }
729
730         
731     //printf("\n\n\n\n");
732     gAlice->TreeR()->Fill();
733   TClonesArray *fRec;
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);
738   }
739   pRICH->ResetRecHits3D();
740
741   free_i3tensor(point,0,kDimensionTheta,0,kDimensionPhi,0,kDimensionOmega);
742   free_i3tensor(point1,0,kDimensionTheta,0,kDimensionPhi,0,kDimensionOmega);
743 }
744
745
746
747 Float_t AliRICHDetect:: Area(Float_t theta,Float_t omega)
748 {
749
750 //
751 // Calculates area of an ellipse for given incidence angles    
752
753
754     Float_t area;
755     const Float_t kHeight=9.25;                       //Distance from Radiator to Pads in pads
756     
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);
758     
759     return (area);
760 }
761
762
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] 
765 {
766   long i,j,nrow=nrh-nrl+1,ncol=nch-ncl+1,ndep=ndh-ndl+1;
767   Int_t ***t;
768   
769   int NR_END=1; 
770
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()");
774   t += NR_END;
775   t -= nrl;
776   
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()");
780   t[nrl] += NR_END;
781   t[nrl] -= ncl;
782   
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;
787   t[nrl][ncl] -= ndl;
788   
789   for(j=ncl+1;j<=nch;j++) t[nrl][j]=t[nrl][j-1]+ndep;
790   for(i=nrl+1;i<=nrh;i++) {
791     t[i]=t[i-1]+ncol;
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;
794   }
795   
796   // return pointer to array of pointers to rows 
797   return t;
798 }
799
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()
802 {
803   int NR_END=1; 
804
805   free((char*) (t[nrl][ncl]+ndl-NR_END));
806   free((char*) (t[nrl]+ncl-NR_END));
807   free((char*) (t+nrl-NR_END));
808 }
809
810
811 Float_t AliRICHDetect:: SnellAngle(Float_t iangle)
812
813
814 // Compute the Snell angle
815
816   Float_t nfreon  = 1.295;
817   Float_t nquartz = 1.585;
818   Float_t ngas = 1;
819
820   Float_t sinrangle;
821   Float_t rangle;
822   Float_t a1, a2;
823
824   a1=nfreon/nquartz;
825   a2=nquartz/ngas;
826
827   sinrangle = a1*a2*sin(iangle);
828
829   if(sinrangle>1.) {
830      rangle = 999.;
831      return rangle;
832   }
833   
834   rangle = asin(sinrangle);  
835   return rangle;
836 }
837
838 Float_t AliRICHDetect:: InvSnellAngle(Float_t rangle)
839
840
841 // Compute the inverse Snell angle
842
843   Float_t nfreon  = 1.295;
844   Float_t nquartz = 1.585;
845   Float_t ngas = 1;
846
847   Float_t siniangle;
848   Float_t iangle;
849   Float_t a1,a2;
850
851   a1=nfreon/nquartz;
852   a2=nquartz/ngas;
853
854   siniangle = sin(rangle)/(a1*a2);
855   iangle = asin(siniangle);
856
857   if(siniangle>1.) {
858      iangle = 999.;
859      return iangle;
860   }
861   
862   iangle = asin(siniangle);
863   return iangle;
864 }
865
866
867
868 //________________________________________________________________________________
869 void AliRICHDetect::CreatePoints(Float_t theta, Float_t phi, Float_t omega, Float_t h)
870 {
871   
872   // Create points along the ellipse equation
873   
874   Int_t s1,s2;
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();
877   
878   static TH2F *REllipse = new TH2F("REllipse","Reconstructed ellipses",150,-25,25,150,-25,25);
879
880   for(Float_t i=0;i<1000;i++)
881     {
882       
883       Float_t counter=0;
884       
885       c0=0;c1=0;c2=0;
886       while((c1*c1-4*c2*c0)<=0 && counter<1000)
887         {
888           //Choose which side to go...
889           if(i>250 && i<750) s1=1; 
890           //if (gRandom->Rndm(1)>.5) s1=1;
891           else s1=-1;
892           //printf("s1:%d\n",s1);
893           //Trial a y
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);
896           Float_t alfa1=theta;
897           Float_t theta1=phi;
898           Float_t omega1=omega;
899           
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);
906           //i+=.01;
907           counter +=1;
908         }
909       
910       if (counter>=1000)
911         y=0; 
912
913       //Choose which side to go...
914       //if(gRandom->Rndm(1)>.5) s=1; 
915       //else s=-1;
916       if(i>500) s2=1;
917       //if (gRandom->Rndm(1)>.5) s2=1;
918       else 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);
922
923       REllipse->Fill(xtrial,y);
924
925       //printf("Coordinates: %f %f %f\n",vectorGlob[0],vectorGlob[1],vectorGlob[2]);
926     }
927   
928   fc3->cd(2);
929   REllipse->Draw();
930 }
931
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)
933 {
934
935   Float_t x,y,y1,h,omega1,omega2;
936
937   Float_t a,b, offset;
938   Float_t a1,b1, offset1;
939
940   h = height;
941
942   //refraction calculation
943   
944   theta = SnellAngle(theta);
945   phi = phi - TMath::Pi();
946   omega1 = SnellAngle(maxOmega);
947   omega2 = SnellAngle(minOmega);
948
949   //maximum zone
950   a = h*(tan(omega1+theta)+tan(omega1-theta))/2;
951   b = h*tan(omega1);
952  
953   offset = h*(tan(omega1+theta)-tan(omega1-theta))/2;
954   
955  
956   //minimum zone
957   a1 = h*(tan(omega2+theta)+tan(omega2-theta))/2;
958   b1 = h*tan(omega2);
959  
960   offset1 = h*(tan(omega2+theta)-tan(omega2-theta))/2;
961   
962
963   //rotation to phi=0
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;
967
968   
969   if(x*x/a+y*y/b<1 && x*x/a1+y1*y1/b1>1)
970     return 1;
971   else
972     return 0;
973
974 }