]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RICH/AliRICH.cxx
309253292b1daf5707535435446ae35ea2d74794
[u/mrichter/AliRoot.git] / RICH / AliRICH.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 #include "AliRICH.h"
17 #include "AliRICHParam.h"
18 #include "AliRICHChamber.h"
19 #include "AliRICHClusterFinder.h"
20 #include <TArrayF.h>
21 #include <TGeometry.h>
22 #include <TBRIK.h>
23 #include <TTUBE.h>
24 #include <TFile.h>
25 #include <TNode.h> 
26 #include <TObjArray.h>
27 #include <TParticle.h>
28 #include <AliStack.h>
29 #include <AliMagF.h>
30 #include <AliRun.h>
31 #include <AliRunDigitizer.h>
32 #include <AliMC.h>
33 #include <TVirtualMC.h>
34  
35 ClassImp(AliRICHhit)
36 //__________________________________________________________________________________________________
37 void AliRICHhit::Print(Option_t*)const
38 {
39   ::Info("hit","Ch=%1i, TID=%6i, eloss=%9.3f eV, in-out dist=%9.4f, OUT(%7.2f,%7.2f,%7.2f)"
40       ,fChamber,fTrack,fEloss*1e9,Length(),fOutX3.X(),fOutX3.Y(),fOutX3.Z());
41 }
42 //__________________________________________________________________________________________________
43 ClassImp(AliRICHdigit)
44 //__________________________________________________________________________________________________
45 void AliRICHdigit::Print(Option_t*)const
46 {
47   ::Info("digit","csxy=%6i, cfm=%9i, c=%2i, x=%3i, y=%3i, q=%8.3f, TID1=%5i, TID2=%5i, TID3=%5i",
48                   Id(),fChFbMip,fChamber,fPadX,fPadY,fQdc,fTracks[0],fTracks[1],fTracks[2]);
49 }
50 //__________________________________________________________________________________________________
51 ClassImp(AliRICHcluster)
52 //__________________________________________________________________________________________________
53 void AliRICHcluster::Print(Option_t*)const
54 {
55   ::Info("cluster","CombiPid=%10i, c=%2i, size=%6i, dim=%5i, x=%7.3f, y=%7.3f, Q=%6i, st=%i",
56            fCombiPid,fChamber,fSize,fDimXY,fX,fY,fQdc,fStatus);
57 }
58 //__________________________________________________________________________________________________
59 ClassImp(AliRICHreco)
60 //__________________________________________________________________________________________________
61 void AliRICHreco::Print(Option_t*)const
62 {
63   ::Info("reco","ThetaCherenkov=%9.6f, Nphotons=%4i, TID=%9i",fThetaCherenkov,fNphotons,fTid);
64 }
65 //__________________________________________________________________________________________________
66 ClassImp(AliRICH)    
67 //__________________________________________________________________________________________________
68 // RICH manager class   
69 //BEGIN_HTML
70 /*
71   <img src="gif/alirich.gif">
72 */
73 //END_HTML
74 //__________________________________________________________________________________________________
75 AliRICH::AliRICH()
76         :AliDetector() 
77 {
78 //Default ctor should not contain any new operators
79   fpParam     =0;
80   fChambers   =0;   
81 //AliDetector ctor deals with Hits and Digits  
82   fSdigits    =0; fNsdigits   =0;
83   fDigitsNew  =0; for(int i=0;i<kNCH;i++) fNdigitsNew[i]  =0;
84   fClusters   =0; for(int i=0;i<kNCH;i++) fNclusters[i]=0;
85   fRecos      =0; fNrecos     =0;
86 }//AliRICH::AliRICH()
87 //__________________________________________________________________________________________________
88 AliRICH::AliRICH(const char *name, const char *title)
89         :AliDetector(name,title)
90 {
91 //Named ctor
92   if(GetDebug())Info("named ctor","Start.");
93   fpParam     =   new AliRICHParam;
94   fChambers = 0;  CreateChambers();
95 //AliDetector ctor deals with Hits and Digits (reset them to 0, does not create them)
96   fHits=       0;     CreateHits();          gAlice->GetMCApp()->AddHitList(fHits);
97   fSdigits=    0;
98   fDigitsNew=  0;
99   fClusters=   0;
100   fRecos      =0;
101   if(GetDebug())Info("named ctor","Stop.");
102 }//AliRICH::AliRICH(const char *name, const char *title)
103 //__________________________________________________________________________________________________
104 AliRICH::~AliRICH()
105 {
106 //dtor
107   if(GetDebug()) Info("dtor","Start.");
108
109   if(fpParam)    delete fpParam;
110   if(fChambers)  delete fChambers;
111   
112   if(fHits)      delete fHits;
113   if(fSdigits)   delete fSdigits;
114   if(fDigits)    delete fDigits;
115   if(fDigitsNew) {fDigitsNew->Delete();   delete fDigitsNew;}
116   if(fClusters)  {fClusters->Delete();    delete fClusters;}
117   if(fRecos)     delete fRecos;
118   if(GetDebug()) Info("dtor","Stop.");    
119 }//AliRICH::~AliRICH()
120 //__________________________________________________________________________________________________
121 void AliRICH::Hits2SDigits()
122 {
123 // Create a list of sdigits corresponding to list of hits. Every hit generates one or more sdigits.
124 //   
125   if(GetDebug()) Info("Hit2SDigits","Start.");
126
127   AliLoader * richLoader = GetLoader();
128   AliRunLoader * runLoader = GetLoader()->GetRunLoader();
129
130   for(Int_t iEventN=0;iEventN<GetLoader()->GetRunLoader()->GetAliRun()->GetEventsPerRun();iEventN++){//events loop
131     runLoader->GetEvent(iEventN);
132   
133     if (!richLoader->TreeH()) richLoader->LoadHits();
134     if (!runLoader->TreeE()) runLoader->LoadHeader(); 
135     if (!runLoader->TreeK()) runLoader->LoadKinematics();//from
136     if (!richLoader->TreeS()) richLoader->MakeTree("S"); MakeBranch("S");//to
137           
138     for(Int_t iPrimN=0;iPrimN<GetLoader()->TreeH()->GetEntries();iPrimN++){//prims loop
139       richLoader->TreeH()->GetEntry(iPrimN);
140       for(Int_t iHitN=0;iHitN<Hits()->GetEntries();iHitN++){//hits loop 
141         AliRICHhit *pHit=(AliRICHhit*)Hits()->At(iHitN);                
142         TVector2 x2 = P()->ShiftToWirePos(C(pHit->C())->Glob2Loc(pHit->OutX3()));                
143         Int_t iTotQdc=P()->TotQdc(x2,pHit->Eloss());
144         if(iTotQdc==0) continue;
145         Int_t iPadXmin,iPadXmax,iPadYmin,iPadYmax;
146         P()->Loc2Area(x2,iPadXmin,iPadYmin,iPadXmax,iPadYmax);//determine affected pads
147         if(GetDebug()) Info("Hits2SDigits","left-down=(%i,%i) right-up=(%i,%i)",iPadXmin,iPadYmin,iPadXmax,iPadYmax);
148         for(Int_t iPadY=iPadYmin;iPadY<=iPadYmax;iPadY++)//affected pads loop
149           for(Int_t iPadX=iPadXmin;iPadX<=iPadXmax;iPadX++){
150             Double_t padQdc=iTotQdc*P()->FracQdc(x2,iPadX,iPadY);
151             if(padQdc>0.1) AddSDigit(pHit->C(),iPadX,iPadY,padQdc,
152               runLoader->Stack()->Particle(pHit->GetTrack())->GetPdgCode(),pHit->GetTrack());
153           }//affected pads loop 
154       }//hits loop
155     }//prims loop
156     richLoader->TreeS()->Fill();
157     richLoader->WriteSDigits("OVERWRITE");
158     ResetSDigits();
159   }//events loop  
160   richLoader->UnloadHits();
161   runLoader->UnloadHeader();
162   runLoader->UnloadKinematics();
163   GetLoader()->UnloadSDigits();  
164   if(GetDebug()) Info("Hit2SDigits","Stop.");
165 }//Hits2SDigits()
166 //__________________________________________________________________________________________________
167 void AliRICH::BuildGeometry() 
168 {
169 //Builds a TNode geometry for event display
170   if(GetDebug())Info("BuildGeometry","Start.");
171   
172   TNode *node, *subnode, *top;
173   top=gAlice->GetGeometry()->GetNode("alice");
174   
175   new TBRIK("S_RICH","S_RICH","void",71.09999,11.5,73.15);
176
177   Float_t wid=P()->SectorSizeX();
178   Float_t len=P()->SectorSizeY();
179   new TBRIK("PHOTO","PHOTO","void",wid/2,0.1,len/2);
180   
181   for(int i=1;i<=kNCH;i++){
182     top->cd();
183     node = new TNode(Form("RICH%i",i),Form("RICH%i",i),"S_RICH",C(i)->X(),C(i)->Y(),C(i)->Z(),C(i)->RotMatrixName());
184     node->SetLineColor(kRed);
185     node->cd();
186     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",wid+P()->DeadZone(),5,len/2+P()->DeadZone()/2,"");
187     subnode->SetLineColor(kGreen);
188     fNodes->Add(subnode);
189     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",0,5,len/2+P()->DeadZone()/2,"");
190     subnode->SetLineColor(kGreen);
191     fNodes->Add(subnode);
192     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",-wid-P()->DeadZone(),5,len/2+P()->DeadZone()/2,"");
193     subnode->SetLineColor(kGreen);
194     fNodes->Add(subnode);
195     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",wid+P()->DeadZone(),5,-len/2-P()->DeadZone()/2,"");
196     subnode->SetLineColor(kGreen);
197     fNodes->Add(subnode);
198     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",0,5,-len/2 -P()->DeadZone()/2,"");
199     subnode->SetLineColor(kGreen);
200     fNodes->Add(subnode);
201     subnode = new TNode("PHOTO1","PHOTO1","PHOTO",-wid-P()->DeadZone(),5,-len/2 - P()->DeadZone()/2,"");
202     subnode->SetLineColor(kGreen);
203     fNodes->Add(subnode);
204     fNodes->Add(node);
205   }  
206   if(GetDebug())Info("BuildGeometry","Stop.");    
207 }//void AliRICH::BuildGeometry()
208
209 //______________________________________________________________________________
210 void AliRICH::CreateMaterials()
211 {
212     //
213     // *** DEFINITION OF AVAILABLE RICH MATERIALS *** 
214   
215 #include "Opticals.h"
216         
217   Float_t a=0,z=0,den=0,radl=0,absl=0;
218   Float_t tmaxfd=-10.0, deemax=-0.2, stemax=-0.1,epsil=0.001, stmin=-0.001; 
219   Int_t   isxfld = gAlice->Field()->Integ();
220   Float_t sxmgmx = gAlice->Field()->Max();
221     
222   AliMaterial( 1, "Air     $",a=14.61,z=7.3, den=0.001205,radl=30420.0,absl=67500);//(Air)
223   AliMedium(1, "DEFAULT MEDIUM AIR$", 1, 0, isxfld, sxmgmx, tmaxfd, stemax, deemax, epsil, stmin);
224   
225   AliMaterial( 6, "HON",      a=12.01,z=6.0, den=0.1,     radl=18.8,   absl=0);    //(C)-equivalent radl
226   AliMedium(2, "HONEYCOMB$", 6, 0, isxfld, sxmgmx, tmaxfd, stemax, deemax, epsil, stmin);
227   
228   AliMaterial(16, "CSI",      a=12.01,z=6.0, den=0.1,     radl=18.8,   absl=0);    //CsI-radl equivalent
229   AliMedium(kCSI, "CSI$", 16, 1, isxfld, sxmgmx,tmaxfd, stemax, deemax, epsil, stmin);
230   
231   AliMaterial(11, "GRI",      a=63.54,z=29.0,den=8.96,    radl=1.43,   absl=0);    //anode grid (Cu) 
232   AliMedium(7, "GRID$", 11, 0, isxfld, sxmgmx, tmaxfd, stemax, deemax, epsil, stmin);
233   
234   AliMaterial(50, "ALUM",     a=26.98,z=13.0,den=2.7,     radl=8.9,    absl=0);    //aluminium sheet (Al)
235   AliMedium(10, "ALUMINUM$", 50, 1, isxfld, sxmgmx, tmaxfd, stemax, deemax, epsil, stmin);
236   
237   AliMaterial(31, "COPPER$",  a=63.54,z=29.0,den=8.96,    radl=1.4,    absl=0);    //(Cu)
238   AliMedium(12, "PCB_COPPER", 31, 0, isxfld, sxmgmx, tmaxfd, stemax, deemax, epsil, stmin);
239   
240   Float_t  aQuartz[2]={28.09,16.0};  Float_t  zQuartz[2]={14.00, 8.0};  Float_t  wmatQuartz[2]={1,2};
241   AliMixture (20, "QUA",aQuartz,zQuartz,den=2.64,-2, wmatQuartz);//Quarz (SiO2) - trasnparent 
242   AliMedium(3, "QUARTZ$", 20, 1, isxfld, sxmgmx, tmaxfd, stemax, deemax, epsil, stmin);
243   
244   AliMixture (21, "QUAO",aQuartz, zQuartz, den=2.64, -2, wmatQuartz);//Quarz (SiO2) - opaque
245   AliMedium(8, "QUARTZO$", 21, 1, isxfld, sxmgmx, tmaxfd, stemax, deemax, epsil, stmin);
246   
247   Float_t  aFreon[2]={12,19};  Float_t  zFreon[2]={6,9};  Float_t wmatFreon[2]={6,14};
248   AliMixture (30, "FRE",aFreon,zFreon,den=1.7,-2,wmatFreon);//Freon (C6F14) 
249   AliMedium(4, "FREON$", 30, 1, isxfld, sxmgmx, tmaxfd, stemax, deemax, epsil, stmin);
250   
251   Float_t aMethane[2]={12.01,1}; Float_t zMethane[2]={6,1}; Float_t wmatMethane[2]={1,4};
252   AliMixture (40, "MET", aMethane, zMethane, den=7.17e-4,-2, wmatMethane);//methane (CH4)     
253   AliMedium(5, "METHANE$", 40, 1, isxfld, sxmgmx, tmaxfd, stemax, deemax, epsil, stmin);
254   
255   AliMixture (41, "METG", aMethane, zMethane, den=7.17e-4, -2, wmatMethane);
256   AliMedium(kGAP, "GAP$", 41, 1, isxfld, sxmgmx,tmaxfd, 0.1, -deemax, epsil, -stmin);
257   
258   Float_t aGlass[5]={12.01, 28.09, 16.,   10.8,  23.};
259   Float_t zGlass[5]={ 6.,   14.,    8.,    5.,   11.};
260   Float_t wGlass[5]={ 0.5,  0.105, 0.355, 0.03,  0.01};
261   AliMixture (32, "GLASS",aGlass, zGlass, den=1.74, 5, wGlass);//Glass 50%C+10.5%Si+35.5%O+3% + 1%
262   AliMedium(11, "GLASS", 32, 0, isxfld, sxmgmx, tmaxfd, stemax, deemax, epsil, stmin);
263             
264   Int_t *idtmed = fIdtmed->GetArray()-999;
265   gMC->SetCerenkov(idtmed[1000], kNbins, aPckov, aAbsCH4,    aQeAll, aIdxCH4);
266   gMC->SetCerenkov(idtmed[1001], kNbins, aPckov, aAbsCH4,    aQeAll, aIdxCH4);
267   gMC->SetCerenkov(idtmed[1002], kNbins, aPckov, aAbsSiO2,   aQeAll, aIdxSiO2);
268   gMC->SetCerenkov(idtmed[1003], kNbins, aPckov, aAbsC6F14,  aQeAll, aIdxC6F14);
269   gMC->SetCerenkov(idtmed[1004], kNbins, aPckov, aAbsCH4,    aQeAll, aIdxCH4);
270   gMC->SetCerenkov(idtmed[1005], kNbins, aPckov, aAbsCsI,    aQeCsI, aIdxCH4);
271   gMC->SetCerenkov(idtmed[1006], kNbins, aPckov, aAbsGrid,   aQeAll, aIdxGrid);
272   gMC->SetCerenkov(idtmed[1007], kNbins, aPckov, aAbsOpSiO2, aQeAll, aIdxOpSiO2);
273   gMC->SetCerenkov(idtmed[1008], kNbins, aPckov, aAbsCH4,    aQeAll, aIdxCH4);
274   gMC->SetCerenkov(idtmed[1009], kNbins, aPckov, aAbsGrid,   aQeAll, aIdxGrid);
275   gMC->SetCerenkov(idtmed[1010], kNbins, aPckov, aAbsOpSiO2, aQeAll, aIdxOpSiO2);
276     
277 }//void AliRICH::CreateMaterials()
278 //__________________________________________________________________________________________________
279 Float_t AliRICH::Fresnel(Float_t ene,Float_t pdoti, Bool_t pola)const
280 {
281
282     //ENE(EV), PDOTI=COS(INC.ANG.), PDOTR=COS(POL.PLANE ROT.ANG.)
283     
284     Float_t en[36] = {5.0,5.1,5.2,5.3,5.4,5.5,5.6,5.7,5.8,5.9,6.0,6.1,6.2,
285                       6.3,6.4,6.5,6.6,6.7,6.8,6.9,7.0,7.1,7.2,7.3,7.4,7.5,7.6,7.7,
286                       7.8,7.9,8.0,8.1,8.2,8.3,8.4,8.5};
287     Float_t csin[36] = {2.14,2.21,2.33,2.48,2.76,2.97,2.99,2.59,2.81,3.05,
288                         2.86,2.53,2.55,2.66,2.79,2.96,3.18,3.05,2.84,2.81,2.38,2.11,
289                         2.01,2.13,2.39,2.73,3.08,3.15,2.95,2.73,2.56,2.41,2.12,1.95,
290                         1.72,1.53};
291     Float_t csik[36] = {0.,0.,0.,0.,0.,0.196,0.408,0.208,0.118,0.49,0.784,0.543,
292                         0.424,0.404,0.371,0.514,0.922,1.102,1.139,1.376,1.461,1.253,0.878,
293                         0.69,0.612,0.649,0.824,1.347,1.571,1.678,1.763,1.857,1.824,1.824,
294                         1.714,1.498};
295     Float_t xe=ene;
296     Int_t  j=Int_t(xe*10)-49;
297     Float_t cn=csin[j]+((csin[j+1]-csin[j])/0.1)*(xe-en[j]);
298     Float_t ck=csik[j]+((csik[j+1]-csik[j])/0.1)*(xe-en[j]);
299
300     //FORMULAE FROM HANDBOOK OF OPTICS, 33.23 OR
301     //W.R. HUNTER, J.O.S.A. 54 (1964),15 , J.O.S.A. 55(1965),1197
302
303     Float_t sinin=TMath::Sqrt(1-pdoti*pdoti);
304     Float_t tanin=sinin/pdoti;
305
306     Float_t c1=cn*cn-ck*ck-sinin*sinin;
307     Float_t c2=4*cn*cn*ck*ck;
308     Float_t aO=TMath::Sqrt(0.5*(TMath::Sqrt(c1*c1+c2)+c1));
309     Float_t b2=0.5*(TMath::Sqrt(c1*c1+c2)-c1);
310     
311     Float_t rs=((aO-pdoti)*(aO-pdoti)+b2)/((aO+pdoti)*(aO+pdoti)+b2);
312     Float_t rp=rs*((aO-sinin*tanin)*(aO-sinin*tanin)+b2)/((aO+sinin*tanin)*(aO+sinin*tanin)+b2);
313     
314
315     //CORRECTION FACTOR FOR SURFACE ROUGHNESS
316     //B.J. STAGG  APPLIED OPTICS, 30(1991),4113
317
318     Float_t sigraf=18.;
319     Float_t lamb=1240/ene;
320     Float_t fresn;
321  
322     Float_t  rO=TMath::Exp(-(4*TMath::Pi()*pdoti*sigraf/lamb)*(4*TMath::Pi()*pdoti*sigraf/lamb));
323
324     if(pola)
325     {
326         Float_t pdotr=0.8;                                 //DEGREE OF POLARIZATION : 1->P , -1->S
327         fresn=0.5*(rp*(1+pdotr)+rs*(1-pdotr));
328     }
329     else
330         fresn=0.5*(rp+rs);
331       
332     fresn = fresn*rO;
333     return(fresn);
334 }//Fresnel()
335 //__________________________________________________________________________________________________
336 Float_t AliRICH::AbsoCH4(Float_t x)const
337 {
338 //Evaluate the absorbtion lenght of CH4
339   Float_t sch4[9] = {.12,.16,.23,.38,.86,2.8,7.9,28.,80.};              //MB X 10^22
340   Float_t em[9] = {8.1,8.158,8.212,8.267,8.322,8.378,8.435,8.493,8.55};
341   const Float_t kLoschmidt=2.686763e19;                                      // LOSCHMIDT NUMBER IN CM-3
342   const Float_t kPressure=750.,kTemperature=283.;                                      
343   const Float_t kPn=kPressure/760.;
344   const Float_t kTn=kTemperature/273.16;
345   const Float_t kC0=-1.655279e-1;
346   const Float_t kC1=6.307392e-2;
347   const Float_t kC2=-8.011441e-3;
348   const Float_t kC3=3.392126e-4;
349                 
350   Float_t crossSection=0;                        
351   if (x<7.75) 
352     crossSection=.06e-22;
353   else if(x>=7.75 && x<=8.1){                 //------ METHANE CROSS SECTION cm-2 ASTROPH. J. 214, L47 (1978)                                               
354         crossSection=(kC0+kC1*x+kC2*x*x+kC3*x*x*x)*1.e-18;
355   }else if (x> 8.1){
356     Int_t j=0;
357     while (x<=em[j] || x>=em[j+1]){
358       j++;
359       Float_t a=(sch4[j+1]-sch4[j])/(em[j+1]-em[j]);
360       crossSection=(sch4[j]+a*(x-em[j]))*1e-22;
361     }
362   }//if
363     
364     Float_t density=kLoschmidt*kPn/kTn; //CH4 molecular density 1/cm-3
365     return 1./(density*crossSection);
366 }//AbsoCH4()
367 //__________________________________________________________________________________________________
368 void AliRICH::MakeBranch(Option_t* option)
369 {
370 //Create Tree branches for the RICH.
371   if(GetDebug())Info("MakeBranch","Start with option= %s.",option);
372     
373   const Int_t kBufferSize = 4000;
374       
375   const char *cH = strstr(option,"H");
376   const char *cD = strstr(option,"D");
377   const char *cR = strstr(option,"R");
378   const char *cS = strstr(option,"S");
379
380   if(cH&&TreeH()){//H
381     CreateHits();      //branch will be created in AliDetector::MakeBranch
382   }//H     
383   AliDetector::MakeBranch(option);//this is after cH because we need to guarantee that fHits array is created
384       
385   if(cS&&fLoader->TreeS()){//S  
386     CreateSDigits();   MakeBranchInTree(fLoader->TreeS(),"RICH",&fSdigits,kBufferSize,0) ;
387   }//S
388    
389   if(cD&&fLoader->TreeD()){//D
390     CreateDigits();
391     for(Int_t i=0;i<kNCH;i++){ 
392       MakeBranchInTree(fLoader->TreeD(),Form("%s%d",GetName(),i+1),&((*fDigitsNew)[i]),kBufferSize,0);
393     }
394   }//D
395   
396   if(cR&&fLoader->TreeR()){//R
397     CreateClusters();
398     for(Int_t i=0;i<kNCH;i++)
399       MakeBranchInTree(fLoader->TreeR(),Form("%sClusters%d",GetName(),i+1), &((*fClusters)[i]), kBufferSize, 0);    
400   }//R
401   if(GetDebug())Info("MakeBranch","Stop.");   
402 }//void AliRICH::MakeBranch(Option_t* option)
403 //__________________________________________________________________________________________________
404 void AliRICH::SetTreeAddress()
405 {
406 //Set branch address for the Hits and Digits Tree.
407   if(GetDebug())Info("SetTreeAddress","Start.");
408       
409   TBranch *branch;
410     
411   if(fLoader->TreeH()){//H
412     if(GetDebug())Info("SetTreeAddress","tree H is requested.");
413     CreateHits();//branch map will be in AliDetector::SetTreeAddress    
414   }//H
415   AliDetector::SetTreeAddress();//this is after TreeH because we need to guarantee that fHits array is created
416
417   if(fLoader->TreeS()){//S
418     if(GetDebug())Info("SetTreeAddress","tree S is requested.");
419     branch=fLoader->TreeS()->GetBranch(GetName());        if(branch){CreateSDigits();   branch->SetAddress(&fSdigits);}
420   }//S
421     
422   if(fLoader->TreeD()){//D    
423     if(GetDebug())Info("SetTreeAddress","tree D is requested.");
424     for(int i=0;i<kNCH;i++){      
425       branch=fLoader->TreeD()->GetBranch(Form("%s%d",GetName(),i+1)); 
426       if(branch){CreateDigits(); branch->SetAddress(&((*fDigitsNew)[i]));}
427     }
428   }//D
429     
430   if(fLoader->TreeR()){//R
431     if(GetDebug())Info("SetTreeAddress","tree R is requested.");
432     for(int i=0;i<kNCH;i++){         
433       branch=fLoader->TreeR()->GetBranch(Form("%sClusters%d" ,GetName(),i+1));
434       if(branch){CreateClusters(); branch->SetAddress(&((*fClusters)[i]));}
435     }
436   }//R
437   if(GetDebug())Info("SetTreeAddress","Stop.");
438 }//void AliRICH::SetTreeAddress()
439 //__________________________________________________________________________________________________
440 void AliRICH::Print(Option_t *option)const
441 {
442 //Debug printout
443   TObject::Print(option);
444   P()->Dump();
445   fChambers->Print(option);  
446 }//void AliRICH::Print(Option_t *option)const
447 //__________________________________________________________________________________________________
448 void AliRICH::CreateGeometry()
449 {
450 //Creates detailed geometry simulation (currently GEANT volumes tree)         
451   if(GetDebug())Info("CreateGeometry","Start.");
452 //Opaque quartz thickness
453   Float_t oquaThickness = .5;
454 //CsI dimensions
455   Float_t pcX=P()->PcSizeX();
456   Float_t pcY=P()->PcSizeY();
457   
458   Int_t *idtmed = fIdtmed->GetArray()-999;
459     
460   Int_t i;
461   Float_t zs;
462   Int_t idrotm[1099];
463   Float_t par[3];
464     
465 //External aluminium box 
466   par[0]=68.8;par[1]=13;par[2]=70.86;  gMC->Gsvolu("RICH", "BOX ", idtmed[1009], par, 3);
467 //Air 
468   par[0]=66.3;   par[1] = 13; par[2] = 68.35;      gMC->Gsvolu("SRIC", "BOX ", idtmed[1000], par, 3); 
469 //Air 2 (cutting the lower part of the box)
470   par[0]=1.25;    par[1] = 3;    par[2] = 70.86;   gMC->Gsvolu("AIR2", "BOX ", idtmed[1000], par, 3);
471 //Air 3 (cutting the lower part of the box)
472   par[0]=66.3;    par[1] = 3;  par[2] = 1.2505;    gMC->Gsvolu("AIR3", "BOX ", idtmed[1000], par, 3);
473 //Honeycomb 
474   par[0]=66.3;par[1]=0.188;  par[2] = 68.35;       gMC->Gsvolu("HONE", "BOX ", idtmed[1001], par, 3);
475 //Aluminium sheet 
476   par[0]=66.3;par[1]=0.025;par[2]=68.35;           gMC->Gsvolu("ALUM", "BOX ", idtmed[1009], par, 3);
477   //par[0] = 66.5; par[1] = .025; par[2] = 63.1;
478 //Quartz 
479   par[0]=P()->QuartzWidth()/2;par[1]=P()->QuartzThickness()/2;par[2]=P()->QuartzLength()/2;
480   gMC->Gsvolu("QUAR", "BOX ", idtmed[1002], par, 3);
481 //Spacers (cylinders) 
482   par[0]=0.;par[1]=.5;par[2]=P()->FreonThickness()/2;  gMC->Gsvolu("SPAC", "TUBE", idtmed[1002], par, 3);    
483 //Feet (freon slabs supports)
484   par[0] = .7;  par[1] = .3;  par[2] = 1.9;        gMC->Gsvolu("FOOT", "BOX", idtmed[1009], par, 3);
485 //Opaque quartz 
486   par[0]=P()->QuartzWidth()/2;par[1]= .2;par[2]=P()->QuartzLength()/2;
487   gMC->Gsvolu("OQUA", "BOX ", idtmed[1007], par, 3);
488 //Frame of opaque quartz
489   par[0]=P()->OuterFreonWidth()/2;par[1]=P()->FreonThickness()/2;par[2]=P()->OuterFreonLength()/2; 
490   gMC->Gsvolu("OQF1", "BOX ", idtmed[1007], par, 3);
491   par[0]=P()->InnerFreonWidth()/2;par[1]=P()->FreonThickness()/2;par[2]=P()->InnerFreonLength()/2; 
492   gMC->Gsvolu("OQF2", "BOX ", idtmed[1007], par, 3);
493 //Freon 
494   par[0]=P()->OuterFreonWidth()/2 - oquaThickness;
495   par[1]=P()->FreonThickness()/2;
496   par[2]=P()->OuterFreonLength()/2 - 2*oquaThickness; 
497   gMC->Gsvolu("FRE1", "BOX ", idtmed[1003], par, 3);
498
499   par[0]=P()->InnerFreonWidth()/2 - oquaThickness;
500   par[1]=P()->FreonThickness()/2;
501   par[2]=P()->InnerFreonLength()/2 - 2*oquaThickness; 
502   gMC->Gsvolu("FRE2", "BOX ", idtmed[1003], par, 3);    
503 //Methane 
504   par[0]=pcX/2;par[1]=P()->GapThickness()/2;par[2]=pcY/2;         gMC->Gsvolu("META","BOX ",idtmed[1004], par, 3);
505 //Methane gap 
506   par[0]=pcX/2;par[1]=P()->ProximityGap()/2;par[2]=pcY/2;gMC->Gsvolu("GAP ","BOX ",(*fIdtmed)[kGAP],par,3);
507 //CsI PC
508   par[0]=pcX/2;par[1]=.25;par[2]=pcY/2;  gMC->Gsvolu("CSI ", "BOX ", (*fIdtmed)[kCSI], par, 3);
509 //Anode grid 
510   par[0] = 0.;par[1] = .001;par[2] = 20.;  gMC->Gsvolu("GRID", "TUBE", idtmed[1006], par, 3);
511
512 //Wire supports
513 //Bar of metal
514   par[0]=pcX/2;par[1]=1.05;par[2]=1.05;  gMC->Gsvolu("WSMe", "BOX ", idtmed[1009], par, 3);
515 //Ceramic pick up (base)
516   par[0]=pcX/2;par[1]= .25;par[2]=1.05;  gMC->Gsvolu("WSG1", "BOX ", idtmed[1010], par, 3);
517 //Ceramic pick up (head)
518   par[0] = pcX/2;par[1] = .1;par[2] = .1;  gMC->Gsvolu("WSG2", "BOX ", idtmed[1010], par, 3);
519
520 //Aluminium supports for methane and CsI
521 //Short bar
522   par[0]=pcX/2;par[1]=P()->GapThickness()/2 + .25; par[2] = (68.35 - pcY/2)/2;
523   gMC->Gsvolu("SMSH", "BOX", idtmed[1009], par, 3);
524 //Long bar
525   par[0]=(66.3 - pcX/2)/2;par[1]=P()->GapThickness()/2+.25;par[2]=pcY/2+68.35-pcY/2;
526   gMC->Gsvolu("SMLG", "BOX", idtmed[1009], par, 3);
527     
528 //Aluminium supports for freon
529 //Short bar
530   par[0] = P()->QuartzWidth()/2; par[1] = .3; par[2] = (68.35 - P()->QuartzLength()/2)/2;
531   gMC->Gsvolu("SFSH", "BOX", idtmed[1009], par, 3);    
532 //Long bar
533   par[0] = (66.3 - P()->QuartzWidth()/2)/2; par[1] = .3;
534   par[2] = P()->QuartzLength()/2 + 68.35 - P()->QuartzLength()/2;
535   gMC->Gsvolu("SFLG", "BOX", idtmed[1009], par, 3);    
536 //PCB backplane
537   par[0] = pcX/2;par[1] = .25; par[2] = pcY/4 -.5025;  gMC->Gsvolu("PCB ", "BOX", idtmed[1011], par, 3);
538
539 //Backplane supports
540 //Aluminium slab
541   par[0] = 33.15;par[1] = 2;par[2] = 21.65;  gMC->Gsvolu("BACK", "BOX", idtmed[1009], par, 3);    
542 //Big hole
543   par[0] = 9.05; par[1] = 2; par[2] = 4.4625;  gMC->Gsvolu("BKHL", "BOX", idtmed[1000], par, 3);
544 //Small hole
545   par[0] = 5.7;par[1] = 2;par[2] = 4.4625;  gMC->Gsvolu("BKHS", "BOX", idtmed[1000], par, 3);
546 //Place holes inside backplane support
547   gMC->Gspos("BKHS", 1, "BACK", .8 + 5.7,0., .6 + 4.4625, 0, "ONLY");
548   gMC->Gspos("BKHS", 2, "BACK", -.8 - 5.7,0., .6 + 4.4625, 0, "ONLY");
549   gMC->Gspos("BKHS", 3, "BACK", .8 + 5.7,0., -.6 - 4.4625, 0, "ONLY");
550   gMC->Gspos("BKHS", 4, "BACK", -.8 - 5.7,0., -.6 - 4.4625, 0, "ONLY");
551   gMC->Gspos("BKHS", 5, "BACK", .8 + 5.7,0., .6 + 8.925 + 1.2 + 4.4625, 0, "ONLY");
552   gMC->Gspos("BKHS", 6, "BACK", -.8 - 5.7,0., .6 + 8.925 + 1.2 + 4.4625, 0, "ONLY");
553   gMC->Gspos("BKHS", 7, "BACK", .8 + 5.7,0., -.6 - 8.925 - 1.2 - 4.4625, 0, "ONLY");
554   gMC->Gspos("BKHS", 8, "BACK", -.8 - 5.7,0., -.6 - 8.925 - 1.2 - 4.4625, 0, "ONLY");
555   gMC->Gspos("BKHL", 1, "BACK", .8 + 11.4 + 1.6 + 9.05, 0., .6 + 4.4625, 0, "ONLY");
556   gMC->Gspos("BKHL", 2, "BACK", -.8 - 11.4 - 1.6 - 9.05, 0., .6 + 4.4625, 0, "ONLY");
557   gMC->Gspos("BKHL", 3, "BACK", .8 + 11.4 + 1.6 + 9.05, 0., -.6 - 4.4625, 0, "ONLY");
558   gMC->Gspos("BKHL", 4, "BACK", -.8 - 11.4 - 1.6 - 9.05, 0., -.6 - 4.4625, 0, "ONLY");
559   gMC->Gspos("BKHL", 5, "BACK", .8 + 11.4+ 1.6 + 9.05, 0., .6 + 8.925 + 1.2 + 4.4625, 0, "ONLY");
560   gMC->Gspos("BKHL", 6, "BACK", -.8 - 11.4 - 1.6 - 9.05, 0., .6 + 8.925 + 1.2 + 4.4625, 0, "ONLY");
561   gMC->Gspos("BKHL", 7, "BACK", .8 + 11.4 + 1.6 + 9.05, 0., -.6 - 8.925 - 1.2 - 4.4625, 0, "ONLY");
562   gMC->Gspos("BKHL", 8, "BACK", -.8 - 11.4 - 1.6 - 9.05, 0., -.6 - 8.925 - 1.2 - 4.4625, 0, "ONLY");
563 //Place material inside RICH 
564   gMC->Gspos("SRIC", 1, "RICH", 0.,0., 0., 0, "ONLY");
565   gMC->Gspos("AIR2", 1, "RICH", 66.3 + 1.2505, 1.276-P()->GapThickness()/2-P()->QuartzThickness()-P()->FreonThickness()- .4 - .6 - .05 - .376 -.5 - 3.35, 0., 0, "ONLY");
566   gMC->Gspos("AIR2", 2, "RICH", -66.3 - 1.2505,1.276-P()->GapThickness()/2-P()->QuartzThickness()-P()->FreonThickness()- .4 - .6 - .05 - .376 -.5 - 3.35, 0., 0, "ONLY");
567   gMC->Gspos("AIR3", 1, "RICH", 0., 1.276-P()->GapThickness()/2 - P()->QuartzThickness() - P()->FreonThickness()- .4 - .6 - .05 - .376 -.5 - 3.35, -68.35 - 1.25, 0, "ONLY");
568   gMC->Gspos("AIR3", 2, "RICH", 0., 1.276 - P()->GapThickness()/2 - P()->QuartzThickness() - P()->FreonThickness()- .4 - .6 - .05 - .376 -.5 - 3.35,  68.35 + 1.25, 0, "ONLY");
569   gMC->Gspos("ALUM", 1, "SRIC", 0., 1.276 - P()->GapThickness()/2 - P()->QuartzThickness() - P()->FreonThickness()- .4 - .6 - .05 - .376 -.025, 0., 0, "ONLY");
570   gMC->Gspos("HONE", 1, "SRIC", 0., 1.276- P()->GapThickness()/2  - P()->QuartzThickness() - P()->FreonThickness()- .4 - .6 - .05 - .188, 0., 0, "ONLY");
571   gMC->Gspos("ALUM", 2, "SRIC", 0., 1.276 - P()->GapThickness()/2 - P()->QuartzThickness() - P()->FreonThickness()- .4 - .6 - .025, 0., 0, "ONLY");
572   gMC->Gspos("FOOT", 1, "SRIC", 64.95, 1.276 - P()->GapThickness()/2 - P()->QuartzThickness() - P()->FreonThickness()- .4 - .3, 36.9, 0, "ONLY");
573   gMC->Gspos("FOOT", 2, "SRIC", 21.65, 1.276 - P()->GapThickness()/2 - P()->QuartzThickness() - P()->FreonThickness()- .4 - .3 , 36.9, 0, "ONLY");
574   gMC->Gspos("FOOT", 3, "SRIC", -21.65, 1.276 - P()->GapThickness()/2 - P()->QuartzThickness() - P()->FreonThickness()- .4 - .3, 36.9, 0, "ONLY");
575   gMC->Gspos("FOOT", 4, "SRIC", -64.95, 1.276 - P()->GapThickness()/2 - P()->QuartzThickness() - P()->FreonThickness()- .4 - .3, 36.9, 0, "ONLY");
576   gMC->Gspos("FOOT", 5, "SRIC", 64.95, 1.276 - P()->GapThickness()/2 - P()->QuartzThickness() - P()->FreonThickness()- .4 - .3, -36.9, 0, "ONLY");
577   gMC->Gspos("FOOT", 6, "SRIC", 21.65, 1.276 - P()->GapThickness()/2 - P()->QuartzThickness() - P()->FreonThickness()- .4 - .3, -36.9, 0, "ONLY");
578   gMC->Gspos("FOOT", 7, "SRIC", -21.65, 1.276 - P()->GapThickness()/2 - P()->QuartzThickness() - P()->FreonThickness()- .4 - .3, -36.9, 0, "ONLY");
579   gMC->Gspos("FOOT", 8, "SRIC", -64.95, 1.276 - P()->GapThickness()/2 - P()->QuartzThickness() - P()->FreonThickness()- .4 - .3, -36.9, 0, "ONLY");
580   gMC->Gspos("OQUA", 1, "SRIC", 0., 1.276 - P()->GapThickness()/2 - P()->QuartzThickness() - P()->FreonThickness()- .2, 0., 0, "ONLY");
581 // Methane supports
582   gMC->Gspos("SMLG", 1, "SRIC", pcX/2 + (66.3 - pcX/2)/2, 1.276 + .25, 0., 0, "ONLY");
583   gMC->Gspos("SMLG", 2, "SRIC", - pcX/2 - (66.3 - pcX/2)/2, 1.276 + .25, 0., 0, "ONLY");
584   gMC->Gspos("SMSH", 1, "SRIC", 0., 1.276 + .25, pcY/2 + (68.35 - pcY/2)/2, 0, "ONLY");
585   gMC->Gspos("SMSH", 2, "SRIC", 0., 1.276 + .25, - pcY/2 - (68.35 - pcY/2)/2, 0, "ONLY");
586 //Freon supports
587   Float_t suppY = 1.276 - P()->GapThickness()/2- P()->QuartzThickness() -P()->FreonThickness() - .2 + .3; //y position of freon supports
588   gMC->Gspos("SFLG", 1, "SRIC", P()->QuartzWidth()/2 + (66.3 - P()->QuartzWidth()/2)/2, suppY, 0., 0, "ONLY");
589   gMC->Gspos("SFLG", 2, "SRIC", - P()->QuartzWidth()/2 - (66.3 - P()->QuartzWidth()/2)/2, suppY, 0., 0, "ONLY");
590   gMC->Gspos("SFSH", 1, "SRIC", 0., suppY, P()->QuartzLength()/2 + (68.35 - P()->QuartzLength()/2)/2, 0, "ONLY");
591   gMC->Gspos("SFSH", 2, "SRIC", 0., suppY, - P()->QuartzLength()/2 - (68.35 - P()->QuartzLength()/2)/2, 0, "ONLY");
592   AliMatrix(idrotm[1019], 0., 0., 90., 0., 90., 90.);
593 //Place spacers
594   Int_t nspacers = 30;
595   for (i = 0; i < nspacers/3; i++) {
596     zs = -11.6/2 + (TMath::Abs(nspacers/6) - i) * 12.2;
597     gMC->Gspos("SPAC", i, "FRE1", 10.5, 0., zs, idrotm[1019], "ONLY");  //Original settings 
598   }
599   for (i = nspacers/3; i < (nspacers*2)/3; i++) {
600     zs = -11.6/2 + (nspacers/3 + TMath::Abs(nspacers/6) - i) * 12.2;
601     gMC->Gspos("SPAC", i, "FRE1", 0, 0., zs, idrotm[1019], "ONLY");  //Original settings 
602   }
603   for (i = (nspacers*2)/3; i < nspacers; ++i) {
604     zs = -11.6/2 + ((nspacers*2)/3 + TMath::Abs(nspacers/6) - i) * 12.2;
605     gMC->Gspos("SPAC", i, "FRE1", -10.5, 0., zs, idrotm[1019], "ONLY"); //Original settings  
606   }
607   for (i = 0; i < nspacers/3; i++) {
608     zs = -11.6/2 + (TMath::Abs(nspacers/6) - i) * 12.2;
609     gMC->Gspos("SPAC", i, "FRE2", 10.5, 0., zs, idrotm[1019], "ONLY");  //Original settings 
610   }
611   for (i = nspacers/3; i < (nspacers*2)/3; i++) {
612     zs = -11.6/2 + (nspacers/3 + TMath::Abs(nspacers/6) - i) * 12.2;
613     gMC->Gspos("SPAC", i, "FRE2", 0, 0., zs, idrotm[1019], "ONLY");  //Original settings 
614   }
615   for (i = (nspacers*2)/3; i < nspacers; ++i) {
616     zs = -11.6/2 + ((nspacers*2)/3 + TMath::Abs(nspacers/6) - i) * 12.2;
617     gMC->Gspos("SPAC", i, "FRE2", -10.5, 0., zs, idrotm[1019], "ONLY"); //Original settings  
618   }
619   gMC->Gspos("FRE1", 1, "OQF1", 0., 0., 0., 0, "ONLY");
620   gMC->Gspos("FRE2", 1, "OQF2", 0., 0., 0., 0, "ONLY");
621   gMC->Gspos("OQF1", 1, "SRIC", P()->OuterFreonWidth()/2 + P()->InnerFreonWidth()/2 + 2, 1.276 - P()->GapThickness()/2- P()->QuartzThickness() -P()->FreonThickness()/2, 0., 0, "ONLY"); //Original settings (31.3)
622   gMC->Gspos("OQF2", 2, "SRIC", 0., 1.276 - P()->GapThickness()/2 - P()->QuartzThickness() - P()->FreonThickness()/2, 0., 0, "ONLY");          //Original settings 
623   gMC->Gspos("OQF1", 3, "SRIC", - (P()->OuterFreonWidth()/2 + P()->InnerFreonWidth()/2) - 2, 1.276 - P()->GapThickness()/2 - P()->QuartzThickness() - P()->FreonThickness()/2, 0., 0, "ONLY");       //Original settings (-31.3)
624   gMC->Gspos("QUAR", 1, "SRIC", 0., 1.276 - P()->GapThickness()/2 - P()->QuartzThickness()/2, 0., 0, "ONLY");
625   gMC->Gspos("GAP ", 1, "META", 0., P()->GapThickness()/2 - P()->ProximityGap()/2 - 0.0001, 0., 0, "ONLY");
626   gMC->Gspos("META", 1, "SRIC", 0., 1.276, 0., 0, "ONLY");
627   gMC->Gspos("CSI ", 1, "SRIC", 0., 1.276 + P()->GapThickness()/2 + .25, 0., 0, "ONLY");
628 //Wire support placing
629   gMC->Gspos("WSG2", 1, "GAP ", 0., P()->ProximityGap()/2 - .1, 0., 0, "ONLY");
630   gMC->Gspos("WSG1", 1, "CSI ", 0., 0., 0., 0, "ONLY");
631   gMC->Gspos("WSMe", 1, "SRIC ", 0., 1.276 + P()->GapThickness()/2 + .5 + 1.05, 0., 0, "ONLY");
632 //Backplane placing
633   gMC->Gspos("BACK", 1, "SRIC ", -33.15, 1.276 + P()->GapThickness()/2 + .5 + 2.1 + 2, 43.3, 0, "ONLY");
634   gMC->Gspos("BACK", 2, "SRIC ", 33.15, 1.276 + P()->GapThickness()/2 + .5 + 2.1 + 2 , 43.3, 0, "ONLY");
635   gMC->Gspos("BACK", 3, "SRIC ", -33.15, 1.276 + P()->GapThickness()/2 + .5 + 2.1 + 2, 0., 0, "ONLY");
636   gMC->Gspos("BACK", 4, "SRIC ", 33.15, 1.276 + P()->GapThickness()/2 + .5 + 2.1 + 2, 0., 0, "ONLY");
637   gMC->Gspos("BACK", 5, "SRIC ", 33.15, 1.276 + P()->GapThickness()/2 + .5 + 2.1 + 2, -43.3, 0, "ONLY");
638   gMC->Gspos("BACK", 6, "SRIC ", -33.15, 1.276 + P()->GapThickness()/2 + .5 + 2.1 + 2, -43.3, 0, "ONLY");
639 //PCB placing
640   gMC->Gspos("PCB ", 1, "SRIC ", 0.,  1.276 + P()->GapThickness()/2 + .5 + 1.05, pcX/4 + .5025 + 2.5, 0, "ONLY");
641   gMC->Gspos("PCB ", 2, "SRIC ", 0.,  1.276 + P()->GapThickness()/2 + .5 + 1.05, -pcX/4 - .5025 - 2.5, 0, "ONLY");
642
643 //place chambers into mother volume ALIC
644   for(int i=1;i<=kNCH;i++){
645     AliMatrix(idrotm[1000+i],C(i)->ThetaXd(),C(i)->PhiXd(),
646                              C(i)->ThetaYd(),C(i)->PhiYd(),
647                              C(i)->ThetaZd(),C(i)->PhiZd());
648     gMC->Gspos("RICH",i,"ALIC",C(i)->X(),C(i)->Y(),C(i)->Z(),idrotm[1000+i], "ONLY");
649   }
650
651   if(GetDebug())Info("CreateGeometry","Stop.");  
652 }//void AliRICH::CreateGeometry()
653 //__________________________________________________________________________________________________
654 void AliRICH::CreateChambers()
655 {
656 //create all RICH Chambers on each call. Previous chambers deleted
657   if(fChambers) delete fChambers;
658   if(GetDebug())Info("CreateChambers","Creating RICH chambers.");
659   fChambers=new TObjArray(kNCH);
660   fChambers->SetOwner();
661   for(int i=0;i<kNCH;i++)  fChambers->AddAt(new AliRICHChamber(i+1,P()),i);  
662 }//void AliRICH::CreateChambers()
663 //__________________________________________________________________________________________________
664 void AliRICH::GenerateFeedbacks(Int_t iChamber,Float_t eloss)
665 {
666 // Generate FeedBack photons 
667 // eloss=0 means photon so only pulse height distribution is to be analysed. This one is done in AliRICHParam::TotQdc()
668   
669   TLorentzVector x4;
670   gMC->TrackPosition(x4);  
671   TVector2 x2=C(iChamber)->Glob2Loc(x4);
672   Int_t sector=P()->Sector(x2);  if(sector==kBad) return; //hit in dead zone nothing to produce
673   Int_t iTotQdc=P()->TotQdc(x2,eloss);
674   Int_t iNphotons=gMC->GetRandom()->Poisson(P()->AlphaFeedback(sector)*iTotQdc);    
675   if(GetDebug())Info("GenerateFeedbacks","N photons=%i",iNphotons);
676   Int_t j;
677   Float_t cthf, phif, enfp = 0, sthf, e1[3], e2[3], e3[3], vmod, uswop,dir[3], phi,pol[3], mom[4];
678 //Generate photons
679   for(Int_t i=0;i<iNphotons;i++){//feedbacks loop
680     Double_t ranf[2];
681     gMC->GetRandom()->RndmArray(2,ranf);    //Sample direction
682     cthf=ranf[0]*2-1.0;
683     if(cthf<0) continue;
684     sthf = TMath::Sqrt((1 - cthf) * (1 + cthf));
685     phif = ranf[1] * 2 * TMath::Pi();
686     
687     if(Double_t randomNumber=gMC->GetRandom()->Rndm()<=0.57)
688       enfp = 7.5e-9;
689     else if(randomNumber<=0.7)
690       enfp = 6.4e-9;
691     else
692       enfp = 7.9e-9;
693     
694
695     dir[0] = sthf * TMath::Sin(phif);    dir[1] = cthf;    dir[2] = sthf * TMath::Cos(phif);
696     gMC->Gdtom(dir, mom, 2);
697     mom[0]*=enfp;    mom[1]*=enfp;    mom[2]*=enfp;
698     mom[3] = TMath::Sqrt(mom[0]*mom[0]+mom[1]*mom[1]+mom[2]*mom[2]);
699     
700     // Polarisation
701     e1[0]=      0;    e1[1]=-dir[2];    e1[2]= dir[1];
702     e2[0]=-dir[1];    e2[1]= dir[0];    e2[2]=      0;
703     e3[0]= dir[1];    e3[1]=      0;    e3[2]=-dir[0];
704     
705     vmod=0;
706     for(j=0;j<3;j++) vmod+=e1[j]*e1[j];
707     if (!vmod) for(j=0;j<3;j++) {
708       uswop=e1[j];
709       e1[j]=e3[j];
710       e3[j]=uswop;
711     }
712     vmod=0;
713     for(j=0;j<3;j++) vmod+=e2[j]*e2[j];
714     if (!vmod) for(j=0;j<3;j++) {
715       uswop=e2[j];
716       e2[j]=e3[j];
717       e3[j]=uswop;
718     }
719     
720     vmod=0;  for(j=0;j<3;j++) vmod+=e1[j]*e1[j];  vmod=TMath::Sqrt(1/vmod);  for(j=0;j<3;j++) e1[j]*=vmod;    
721     vmod=0;  for(j=0;j<3;j++) vmod+=e2[j]*e2[j];  vmod=TMath::Sqrt(1/vmod);  for(j=0;j<3;j++) e2[j]*=vmod;
722     
723     phi = gMC->GetRandom()->Rndm()* 2 * TMath::Pi();
724     for(j=0;j<3;j++) pol[j]=e1[j]*TMath::Sin(phi)+e2[j]*TMath::Cos(phi);
725     gMC->Gdtom(pol, pol, 2);
726     Int_t outputNtracksStored;    
727     gAlice->GetMCApp()->PushTrack(1,                 //do not transport
728                      gAlice->GetMCApp()->GetCurrentTrackNumber(),//parent track 
729                      kFeedback,                      //PID
730                      mom[0],mom[1],mom[2],mom[3],    //track momentum  
731                      x4.X(),x4.Y(),x4.Z(),x4.T(),    //track origin 
732                      pol[0],pol[1],pol[2],           //polarization
733                      kPFeedBackPhoton,
734                      outputNtracksStored,
735                      1.0);    
736   }//feedbacks loop
737 }//GenerateFeedbacks()
738 //__________________________________________________________________________________________________
739
740 void AliRICH::Reconstruct() const
741 {
742 // reconstruct clusters
743
744   AliRICHClusterFinder clusterer(const_cast<AliRICH*>(this));
745   clusterer.Exec();
746 }
747