]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDv1.cxx
Changes for correct initialization of Geant4 (Mihaela)
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDv1.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 #include "AliHMPIDv1.h"     //class header
18 #include "AliHMPIDParam.h"  //CreateMaterials()
19 #include "AliHMPIDHit.h"    //Hits2SDigs(),StepManager()
20 #include "AliHMPIDDigit.h"  //CreateMaterials()
21 #include "AliRawReader.h"  //Raw2SDigits()
22 #include <TParticle.h>     //Hits2SDigits()
23 #include <TRandom.h> 
24 #include <TVirtualMC.h>    //StepManager() for gMC
25 #include <TPDGCode.h>      //StepHistory() 
26 #include <AliStack.h>      //StepManager(),Hits2SDigits()
27 #include <AliLoader.h>        //Hits2SDigits()
28 #include <AliRunLoader.h>     //Hits2SDigits()
29 #include <AliConst.h>
30 #include <AliPDG.h>
31 #include <AliMC.h>            //StepManager()      
32 #include <AliRawDataHeader.h> //Digits2Raw()
33 #include <AliDAQ.h>           //Digits2Raw()
34 #include <AliRun.h>           //CreateMaterials()    
35 #include <AliMagF.h>          //CreateMaterials()
36 #include <TGeoManager.h>      //CreateGeometry()
37 #include <TMultiGraph.h>      //Optics() 
38 #include <TGraph.h>           //Optics() 
39 #include <TLegend.h>          //Optics() 
40 #include <TCanvas.h>          //Optics() 
41 #include <TF2.h>              //CreateMaterials()
42 #include <AliCDBManager.h>    //CreateMaterials()
43 #include <AliCDBEntry.h>      //CreateMaterials()
44  
45 ClassImp(AliHMPIDv1)    
46 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
47 void AliHMPIDv1::AddAlignableVolumes()const
48 {
49 // Associates the symbolic volume name with the corresponding volume path. Interface methode from AliModule ivoked from AliMC
50 // Arguments: none
51 //   Returns: none   
52   for(Int_t i=0;i<7;i++)
53     gGeoManager->SetAlignableEntry(Form("/HMPID/Chamber%i",i),Form("ALIC_1/HMPID_%i",i));
54 }
55 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
56 void AliHMPIDv1::CreateMaterials()
57 {
58 // Definition of available HMPID materials  
59 // Arguments: none
60 //   Returns: none    
61   AliDebug(1,"Start v1 HMPID.");
62     
63 //data from PDG booklet 2002     density [gr/cm^3] rad len [cm] abs len [cm]    
64   Float_t   aAir[4]={12,14,16,36}    ,   zAir[4]={6,7,8,18} ,   wAir[4]={0.000124,0.755267,0.231781,0.012827} , dAir=0.00120479; Int_t nAir=4;//mixture 0.9999999
65   Float_t aC6F14[2]={ 12.01 , 18.99} , zC6F14[2]={ 6 , 9}   , wC6F14[2]={6 , 14} , dC6F14=1.68    ; Int_t nC6F14=-2;
66   Float_t  aSiO2[2]={ 28.09 , 15.99} ,  zSiO2[2]={14 , 8}   ,  wSiO2[2]={1 ,  2} ,  dSiO2=2.64    ; Int_t  nSiO2=-2; 
67   Float_t   aCH4[2]={ 12.01 ,  1.01} ,   zCH4[2]={ 6 , 1}   ,   wCH4[2]={1 ,  4} ,   dCH4=7.17e-4 ; Int_t   nCH4=-2; 
68   Float_t   aCsI[2]={132.90 ,126.90} ,   zCsI[2]={55 ,53}   ,   wCsI[2]={1 ,  1} ,   dCsI=0.1     ; Int_t   nCsI=-2; 
69   Float_t     aRoha= 12.01 ,               zRoha=  6 ,                              dRoha=  0.10 , radRoha= 18.80 , absRoha=  86.3/dRoha; //special material- quazi carbon
70   Float_t       aCu= 63.55 ,                 zCu= 29 ,                                dCu=  8.96 ,   radCu=  1.43 ,   absCu= 134.9/dCu  ;
71   Float_t        aW=183.84 ,                  zW= 74 ,                                 dW= 19.30 ,    radW=  0.35 ,    absW= 185.0/dW   ;
72   Float_t       aAl= 26.98 ,                 zAl= 13 ,                                dAl=  2.70 ,   radAl=  8.90 ,   absAl= 106.4/dAl  ;
73   
74     Int_t   matId=0;                           //tmp material id number
75     Int_t   unsens =  0, sens=1;               //sensitive or unsensitive medium
76     Int_t   itgfld = gAlice->Field()->Integ(); //type of field intergration 0 no field -1 user in guswim 1 Runge Kutta 2 helix 3 const field along z
77     Float_t maxfld = gAlice->Field()->Max();   //max field value
78     Float_t tmaxfd = -10.0;                    //max deflection angle due to magnetic field in one step
79     Float_t deemax = - 0.2;                    //max fractional energy loss in one step   
80     Float_t stemax = - 0.1;                    //mas step allowed [cm]
81     Float_t epsil  =   0.001;                  //abs tracking precision [cm]   
82     Float_t stmin  = - 0.001;                  //min step size [cm] in continius process transport, negative value: choose it automatically
83     AliMixture(++matId,"Air"  ,aAir  ,zAir  ,dAir  ,nAir  ,wAir  ); AliMedium(kAir  ,"Air"  ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
84     AliMixture(++matId,"C6F14",aC6F14,zC6F14,dC6F14,nC6F14,wC6F14); AliMedium(kC6F14,"C6F14",matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);      
85     AliMixture(++matId,"SiO2" ,aSiO2 ,zSiO2 ,dSiO2 ,nSiO2 ,wSiO2 ); AliMedium(kSiO2 ,"SiO2" ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);    
86     AliMixture(++matId,"CH4"  ,aCH4  ,zCH4  ,dCH4  ,nCH4  ,wCH4  ); AliMedium(kCH4  ,"CH4"  ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);  
87     AliMixture(++matId,"CsI"  ,aCsI  ,zCsI  ,dCsI  ,nCsI  ,wCsI  ); AliMedium(kCsI  ,"CsI"  ,matId,   sens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);//sensitive
88   
89     AliMaterial(++matId,"Roha",aRoha,zRoha,dRoha,radRoha,absRoha);  AliMedium(kRoha,"Roha", matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
90     AliMaterial(++matId,"Cu"  ,aCu  ,zCu  ,dCu  ,radCu  ,absCu  );  AliMedium(kCu  ,"Cu"  , matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
91     AliMaterial(++matId,"W"   ,aW   ,zW   ,dW   ,radW   ,absW   );  AliMedium(kW   ,"W"   , matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
92     AliMaterial(++matId,"Al"  ,aAl  ,zAl  ,dAl  ,radAl  ,absAl  );  AliMedium(kAl  ,"Al"  , matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
93     
94     AliDebug(1,"Stop v1 HMPID.");
95
96   TString ttl=GetTitle(); if(!ttl.Contains("ShowOptics")) return;              //user didn't aks to plot optical curves
97   
98   const Double_t kWidth=0.25,kHeight=0.2;  
99   const Int_t kRadM=24 , kRadC=kRed;  
100   const Int_t kWinM=26 , kWinC=kBlue;  
101   const Int_t kGapM=25 , kGapC=kGreen;  
102   const Int_t kPcM = 2 , kPcC =kMagenta;  
103   const Int_t kNbins=30;       //number of photon energy points
104
105   Float_t aEckov [kNbins]; 
106   Float_t aAbsRad[kNbins], aAbsWin[kNbins], aAbsGap[kNbins];
107   Float_t aIdxRad[kNbins], aIdxWin[kNbins], aIdxGap[kNbins]; 
108   Float_t aQePc [kNbins];  
109   Float_t aTraRad[kNbins],aTraWin[kNbins],aTraGap[kNbins],aTraTot[kNbins];
110   for(Int_t i=0;i<kNbins;i++){//calculate probability for photon to survive during transversing a volume of material with absorption length  
111     aTraRad[i]=TMath::Exp(-AliHMPIDDigit::SizeRad()/ (aAbsRad[i]+0.0001)); //radiator
112     aTraWin[i]=TMath::Exp(-AliHMPIDDigit::SizeWin()/ (aAbsWin[i] +0.0001)); //window
113     aTraGap[i]=TMath::Exp(-AliHMPIDDigit::SizeGap()/ (aAbsGap[i]  +0.0001)); //from window to PC   
114     aTraTot[i]=aTraRad[i]*aTraWin[i]*aTraGap[i]*aQePc[i];
115   }
116   
117   TGraph *pRaAG=new TGraph(kNbins,aEckov,aAbsRad);pRaAG->SetMarkerStyle(kRadM);pRaAG->SetMarkerColor(kRadC);
118   TGraph *pRaIG=new TGraph(kNbins,aEckov,aIdxRad);pRaIG->SetMarkerStyle(kRadM);pRaIG->SetMarkerColor(kRadC);
119   TGraph *pRaTG=new TGraph(kNbins,aEckov,aTraRad);pRaTG->SetMarkerStyle(kRadM);pRaTG->SetMarkerColor(kRadC);  
120   
121   TGraph *pWiAG=new TGraph(kNbins,aEckov,aAbsWin);pWiAG->SetMarkerStyle(kWinM);pWiAG->SetMarkerColor(kWinC);
122   TGraph *pWiIG=new TGraph(kNbins,aEckov,aIdxWin);pWiIG->SetMarkerStyle(kWinM);pWiIG->SetMarkerColor(kWinC);  
123   TGraph *pWiTG=new TGraph(kNbins,aEckov,aTraWin);pWiTG->SetMarkerStyle(kWinM);pWiTG->SetMarkerColor(kWinC);  
124   
125   TGraph *pGaAG=new TGraph(kNbins,aEckov,aAbsGap);pGaAG->SetMarkerStyle(kGapM);pGaAG->SetMarkerColor(kGapC);
126   TGraph *pGaIG=new TGraph(kNbins,aEckov,aIdxGap);pGaIG->SetMarkerStyle(kGapM);pGaIG->SetMarkerColor(kGapC);
127   TGraph *pGaTG=new TGraph(kNbins,aEckov,aTraGap);pGaTG->SetMarkerStyle(kGapM);pGaTG->SetMarkerColor(kGapC);   
128   
129   TGraph *pQeG =new TGraph(kNbins,aEckov,aQePc);  pQeG  ->SetMarkerStyle(kPcM );pQeG->SetMarkerColor(kPcC);
130   TGraph *pToG =new TGraph(kNbins,aEckov,aTraTot);pToG  ->SetMarkerStyle(30)   ;pToG->SetMarkerColor(kYellow);  
131   
132   TMultiGraph *pIdxMG=new TMultiGraph("idx","Ref index;E_{#check{C}} [GeV]");       
133   TMultiGraph *pAbsMG=new TMultiGraph("abs","Absorption [cm];E_{#check{C}} [GeV]"); 
134   TMultiGraph *pTraMG=new TMultiGraph("tra","Transmission;E_{#check{C}} [GeV]");    TLegend *pTraLe=new TLegend(0.2,0.4,0.2+kWidth,0.4+kHeight);
135   pAbsMG->Add(pRaAG);  pIdxMG->Add(pRaIG);     pTraMG->Add(pRaTG);     pTraLe->AddEntry(pRaTG, "Rad", "p");           
136   pAbsMG->Add(pWiAG);  pIdxMG->Add(pWiIG);     pTraMG->Add(pWiTG);     pTraLe->AddEntry(pWiTG, "Win", "p");               
137   pAbsMG->Add(pGaAG);  pIdxMG->Add(pGaIG);     pTraMG->Add(pGaTG);     pTraLe->AddEntry(pGaTG, "Gap", "p");               
138                                                pTraMG->Add(pToG);      pTraLe->AddEntry(pToG,  "Tot", "p");          
139                                                pTraMG->Add(pQeG);      pTraLe->AddEntry(pQeG,   "QE" , "p");  
140   TCanvas *pC=new TCanvas("c1","HMPID optics to check",1100,900);  pC->Divide(2,2);           
141   pC->cd(1);                    pIdxMG->Draw("AP"); 
142   pC->cd(2);  gPad->SetLogy();  pAbsMG->Draw("AP"); 
143   pC->cd(3);                    pTraLe->Draw();    
144   pC->cd(4);                    pTraMG->Draw("AP");
145 }//void AliHMPID::CreateMaterials()
146 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
147 void AliHMPIDv1::CreateGeometry()
148 {
149 //Creates detailed geometry simulation (currently GEANT volumes tree)         
150   AliDebug(1,"Start main.");
151   if(!gMC->IsRootGeometrySupported()) return;     
152   
153   Double_t cm=1,mm=0.1*cm,mkm=0.001*mm,dx,dy,dz;//default is cm
154   
155   TGeoVolume *pRich=gGeoManager->MakeBox("HMPID",gGeoManager->GetMedium("HMPID_CH4"),dx=(6*mm+1681*mm+6*mm)/2,  //main HMPID volume
156                                                                                    dy=(6*mm+1466*mm+6*mm)/2,
157                                                                                    dz=(80*mm+40*mm)*2/2);     //x,y taken from 2033P1  z from p84 TDR  
158   const Double_t kAngHor=19.5; //  horizontal angle between chambers  19.5 grad
159   const Double_t kAngVer=20;   //  vertical angle between chambers    20   grad     
160   const Double_t kAngCom=30;   //  common HMPID rotation with respect to x axis  30   grad     
161   const Double_t trans[3]={490,0,0}; //center of the chamber is on window-gap surface
162   for(Int_t iCh=0;iCh<7;iCh++){//place 7 chambers
163     TGeoHMatrix *pMatrix=new TGeoHMatrix;
164     pMatrix->RotateY(90);           //rotate around y since initial position is in XY plane -> now in YZ plane
165     pMatrix->SetTranslation(trans); //now plane in YZ is shifted along x 
166     switch(iCh){
167       case 0:                pMatrix->RotateY(kAngHor);  pMatrix->RotateZ(-kAngVer);  break; //right and down 
168       case 1:                                            pMatrix->RotateZ(-kAngVer);  break; //down              
169       case 2:                pMatrix->RotateY(kAngHor);                               break; //right 
170       case 3:                                                                         break; //no rotation
171       case 4:                pMatrix->RotateY(-kAngHor);                              break; //left   
172       case 5:                                            pMatrix->RotateZ(kAngVer);   break; //up
173       case 6:                pMatrix->RotateY(-kAngHor); pMatrix->RotateZ(kAngVer);   break; //left and up 
174     }
175     pMatrix->RotateZ(kAngCom);     //apply common rotation  in XY plane    
176     gGeoManager->GetVolume("ALIC")->AddNode(pRich,iCh,pMatrix);
177   }
178
179   Float_t par[3];
180   Int_t matrixIdReturn=0; //matrix id returned by AliMatrix
181 //Pad Panel frame  6 sectors
182   par[0]=648*mm/2;par[1]=  411*mm/2;par[2]=40  *mm/2;gMC->Gsvolu("Rppf"     ,"BOX ",(*fIdtmed)[kAl]  ,par,3);//PPF 2001P2 inner size of the slab by 1mm more
183   par[0]=181*mm/2;par[1]=89.25*mm/2;par[2]=38.3*mm/2;gMC->Gsvolu("RppfLarge","BOX ",(*fIdtmed)[kAir] ,par,3);//large whole
184   par[0]=114*mm/2;par[1]=89.25*mm/2;par[2]=38.3*mm/2;gMC->Gsvolu("RppfSmall","BOX ",(*fIdtmed)[kAir] ,par,3);//small whole
185   par[0]=644*mm/2;par[1]=  407*mm/2;par[2]= 1.7*mm/2;gMC->Gsvolu("Rpc"      ,"BOX ",(*fIdtmed)[kCsI] ,par,3);//by 0.2 mm more then actual size (PCB 2006P1)
186   
187   gMC->Gspos("Rppf",0,"HMPID",    -335*mm,      -433*mm,  8*cm+20*mm,  0,"ONLY");//F1 2040P1 z p.84 TDR
188   gMC->Gspos("Rppf",1,"HMPID",    +335*mm,      -433*mm,  8*cm+20*mm,  0,"ONLY");
189   gMC->Gspos("Rppf",2,"HMPID",    -335*mm,         0*mm,  8*cm+20*mm,  0,"ONLY");
190   gMC->Gspos("Rppf",3,"HMPID",    +335*mm,         0*mm,  8*cm+20*mm,  0,"ONLY");
191   gMC->Gspos("Rppf",4,"HMPID",    -335*mm,      +433*mm,  8*cm+20*mm,  0,"ONLY");
192   gMC->Gspos("Rppf",5,"HMPID",    +335*mm,      +433*mm,  8*cm+20*mm,  0,"ONLY");  
193     gMC->Gspos("Rpc"      ,1,"Rppf",       0*mm,         0*mm,   -19.15*mm,  0,"ONLY");//PPF 2001P2 
194     gMC->Gspos("RppfLarge",1,"Rppf",  -224.5*mm,  -151.875*mm,     0.85*mm,  0,"ONLY");
195     gMC->Gspos("RppfLarge",2,"Rppf",  -224.5*mm,  - 50.625*mm,     0.85*mm,  0,"ONLY");
196     gMC->Gspos("RppfLarge",3,"Rppf",  -224.5*mm,  + 50.625*mm,     0.85*mm,  0,"ONLY");
197     gMC->Gspos("RppfLarge",4,"Rppf",  -224.5*mm,  +151.875*mm,     0.85*mm,  0,"ONLY");
198     gMC->Gspos("RppfSmall",1,"Rppf",  - 65.0*mm,  -151.875*mm,     0.85*mm,  0,"ONLY");
199     gMC->Gspos("RppfSmall",2,"Rppf",  - 65.0*mm,  - 50.625*mm,     0.85*mm,  0,"ONLY");
200     gMC->Gspos("RppfSmall",3,"Rppf",  - 65.0*mm,  + 50.625*mm,     0.85*mm,  0,"ONLY");
201     gMC->Gspos("RppfSmall",4,"Rppf",  - 65.0*mm,  +151.875*mm,     0.85*mm,  0,"ONLY");
202     gMC->Gspos("RppfSmall",5,"Rppf",  + 65.0*mm,  -151.875*mm,     0.85*mm,  0,"ONLY");
203     gMC->Gspos("RppfSmall",6,"Rppf",  + 65.0*mm,  - 50.625*mm,     0.85*mm,  0,"ONLY");
204     gMC->Gspos("RppfSmall",7,"Rppf",  + 65.0*mm,  + 50.625*mm,     0.85*mm,  0,"ONLY");
205     gMC->Gspos("RppfSmall",8,"Rppf",  + 65.0*mm,  +151.875*mm,     0.85*mm,  0,"ONLY"); 
206     gMC->Gspos("RppfLarge",5,"Rppf",  +224.5*mm,  -151.875*mm,     0.85*mm,  0,"ONLY");
207     gMC->Gspos("RppfLarge",6,"Rppf",  +224.5*mm,  - 50.625*mm,     0.85*mm,  0,"ONLY");
208     gMC->Gspos("RppfLarge",7,"Rppf",  +224.5*mm,  + 50.625*mm,     0.85*mm,  0,"ONLY");
209     gMC->Gspos("RppfLarge",8,"Rppf",  +224.5*mm,  +151.875*mm,     0.85*mm,  0,"ONLY");
210 //Gap - anod wires 6 copies to HMPID
211   par[0]=648*mm/2;par[1]=  411*mm/2 ;par[2]=4.45*mm/2;gMC->Gsvolu("Rgap","BOX ",(*fIdtmed)[kCH4] ,par,3);//xy as PPF 2001P2 z WP 2099P1
212   par[0]=  0*mm  ;par[1]=  20*mkm/2 ;par[2]= 648*mm/2;gMC->Gsvolu("Rano","TUBE",(*fIdtmed)[kW]   ,par,3);//WP 2099P1 z = gap x PPF 2001P2
213   AliMatrix(matrixIdReturn,180,0, 90,90, 90,0); //wires along x
214   
215   gMC->Gspos("Rgap",0,"HMPID",    -335*mm,      -433*mm,8*cm-2.225*mm, 0,"ONLY"); //F1 2040P1 z WP 2099P1
216   gMC->Gspos("Rgap",1,"HMPID",    +335*mm,      -433*mm,8*cm-2.225*mm, 0,"ONLY"); 
217   gMC->Gspos("Rgap",2,"HMPID",    -335*mm,         0*mm,8*cm-2.225*mm, 0,"ONLY"); 
218   gMC->Gspos("Rgap",3,"HMPID",    +335*mm,         0*mm,8*cm-2.225*mm, 0,"ONLY"); 
219   gMC->Gspos("Rgap",4,"HMPID",    -335*mm,      +433*mm,8*cm-2.225*mm, 0,"ONLY"); 
220   gMC->Gspos("Rgap",5,"HMPID",    +335*mm,      +433*mm,8*cm-2.225*mm, 0,"ONLY"); 
221   for(int i=1;i<=96;i++)
222     gMC->Gspos("Rano",i,"Rgap",     0*mm, -411/2*mm+i*4*mm, 0.185*mm, matrixIdReturn,"ONLY"); //WP 2099P1  
223 //Defines radiators geometry  
224   par[0]=1330*mm/2 ;par[1]= 413*mm/2  ;par[2]=  24*mm/2;  gMC->Gsvolu("Rrad"      ,"BOX ",(*fIdtmed)[kC6F14]     ,par,3); // Rad 2011P1
225   par[0]=1330*mm/2 ;par[1]= 413*mm/2  ;par[2]=   4*mm/2;  gMC->Gsvolu("RradFront" ,"BOX ",(*fIdtmed)[kRoha]      ,par,3); //front 
226   par[0]=1330*mm/2 ;par[1]= 413*mm/2  ;par[2]=   5*mm/2;  gMC->Gsvolu("RradWin"   ,"BOX ",(*fIdtmed)[kSiO2]      ,par,3); //window
227   par[0]=1330*mm/2 ;par[1]=   5*mm/2  ;par[2]=  15*mm/2;  gMC->Gsvolu("RradLong"  ,"BOX ",(*fIdtmed)[kRoha]      ,par,3); //long side  
228   par[0]=  10*mm/2 ;par[1]= 403*mm/2  ;par[2]=  15*mm/2;  gMC->Gsvolu("RradShort" ,"BOX ",(*fIdtmed)[kRoha]      ,par,3); //short side 
229   par[0]=   0      ;par[1]=  10*mm/2  ;par[2]=  15*mm/2;  gMC->Gsvolu("RradSpacer","TUBE",(*fIdtmed)[kSiO2]      ,par,3); //spacer        
230     
231   gMC->Gspos("Rrad",1,"HMPID",   0*mm,-434*mm,   -12*mm,  0,"ONLY"); //3 radiators to HMPID
232   gMC->Gspos("Rrad",2,"HMPID",   0*mm,   0*mm,   -12*mm,  0,"ONLY"); 
233   gMC->Gspos("Rrad",3,"HMPID",   0*mm,+434*mm,   -12*mm,  0,"ONLY"); 
234     gMC->Gspos("RradFront",1,"Rrad",   0*mm,   0*mm, -10.0*mm,  0,"ONLY"); //front cover 
235     gMC->Gspos("RradWin"  ,1,"Rrad",   0*mm,   0*mm,   9.5*mm,  0,"ONLY"); //quartz window (back cover)
236     gMC->Gspos("RradLong" ,1,"Rrad",   0*mm,-204*mm,  -0.5*mm,  0,"ONLY"); //long side
237     gMC->Gspos("RradLong" ,2,"Rrad",   0*mm,+204*mm,  -0.5*mm,  0,"ONLY"); //long side
238     gMC->Gspos("RradShort",1,"Rrad",-660*mm,   0*mm,  -0.5*mm,  0,"ONLY"); //short side
239     gMC->Gspos("RradShort",2,"Rrad",+660*mm,   0*mm,  -0.5*mm,  0,"ONLY"); //short side 
240     for(int i=0;i<3;i++)
241       for(int j=0;j<10;j++)
242         gMC->Gspos("RradSpacer",10*i+j,"Rrad",-1330*mm/2+116*mm+j*122*mm,(i-1)*105*mm,-0.5*mm,0,"ONLY");//spacers
243 //Defines SandBox geometry
244   par[0]=1419*mm/2 ;par[1]=1378*mm/2;par[2]=50.5*mm/2; gMC->Gsvolu("Rsb"     ,"BOX ",(*fIdtmed)[kAir]  ,par,3);  //2072P1   
245   par[0]=1419*mm/2 ;par[1]=1378*mm/2;par[2]= 0.5*mm/2; gMC->Gsvolu("RsbCover","BOX ",(*fIdtmed)[kAl]   ,par,3);  //cover
246   par[0]=1359*mm/2 ;par[1]=1318*mm/2;par[2]=49.5*mm/2; gMC->Gsvolu("RsbComb" ,"BOX ",(*fIdtmed)[kRoha] ,par,3);  //honeycomb structure 
247   
248   gMC->Gspos("Rsb",1,"HMPID",   0*mm, 0*mm, -73.75*mm, 0,"ONLY"); //p.84 TDR sandbox to rich
249     gMC->Gspos("RsbComb" ,1,"Rsb", 0*mm, 0*mm,      0*mm, 0,"ONLY"); //2072P1 honeycomv to sandbox
250     gMC->Gspos("RsbCover",1,"Rsb", 0*mm, 0*mm,    +25*mm, 0,"ONLY"); //cover to sandbox
251     gMC->Gspos("RsbCover",2,"Rsb", 0*mm, 0*mm,    -25*mm, 0,"ONLY"); //cover to sandbox
252   AliDebug(1,"Stop v1. HMPID option");  
253 }//CreateGeometry()
254 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
255 void AliHMPIDv1::Init()
256 {
257 // This methode defines ID for sensitive volumes, i.e. such geometry volumes for which there are if(gMC->CurrentVolID()==XXX) statements in StepManager()
258 // Arguments: none
259 //   Returns: none      
260   AliDebug(1,"Start v1 HMPID.");    
261   fIdRad     = gMC->VolId("Rrad");
262   fIdWin     = gMC->VolId("RradWin");
263   fIdPc      = gMC->VolId("Rpc");
264   fIdAmpGap  = gMC->VolId("Rgap");
265   fIdProxGap = gMC->VolId("Rgap");
266   const Int_t kNbins=30;       //number of photon energy points
267   Float_t emin=5.5,emax=8.5;         //Photon energy range,[eV]
268   Float_t aEckov [kNbins]; 
269   Float_t aAbsRad[kNbins], aAbsWin[kNbins], aAbsGap[kNbins], aAbsMet[kNbins];
270   Float_t aIdxRad[kNbins], aIdxWin[kNbins], aIdxGap[kNbins], aIdxMet[kNbins], aIdxPc[kNbins]; 
271   Float_t                                                    aQeAll [kNbins], aQePc [kNbins];
272
273   TF2 *pRaIF=new TF2("RidxRad","sqrt(1+0.554*(1239.84/x)^2/((1239.84/x)^2-5796)-0.0005*(y-20))"                                       ,emin,emax,0,50); //DiMauro mail temp 0-50 degrees C
274   TF1 *pWiIF=new TF1("RidxWin","sqrt(1+46.411/(10.666*10.666-x*x)+228.71/(18.125*18.125-x*x))"                                        ,emin,emax);      //SiO2 idx TDR p.35
275   TF1 *pGaIF=new TF1("RidxGap","1+0.12489e-6/(2.62e-4 - x*x/1239.84/1239.84)"                                                         ,emin,emax);      //?????? from where  
276
277   TF1 *pRaAF=new TF1("RabsRad","(x<7.8)*(gaus+gaus(3))+(x>=7.8)*0.0001"                                                               ,emin,emax);  //fit from DiMauro data 28.10.03 
278   pRaAF->SetParameters(3.20491e16,-0.00917890,0.742402,3035.37,4.81171,0.626309);
279   TF1 *pWiAF=new TF1("RabsWin","(x<8.2)*(818.8638-301.0436*x+36.89642*x*x-1.507555*x*x*x)+(x>=8.2)*0.0001"                            ,emin,emax);  //fit from DiMauro data 28.10.03 
280   TF1 *pGaAF=new TF1("RabsGap","(x<7.75)*6512.399+(x>=7.75)*3.90743e-2/(-1.655279e-1+6.307392e-2*x-8.011441e-3*x*x+3.392126e-4*x*x*x)",emin,emax);  //????? from where  
281   
282   TF1 *pQeF =new TF1("Qe"    ,"0+(x>6.07267)*0.344811*(1-exp(-1.29730*(x-6.07267)))"                                                  ,emin,emax);  //fit from DiMauro data 28.10.03  
283                             
284   for(Int_t i=0;i<kNbins;i++){
285     Float_t eV=emin+0.1*i;  //Ckov energy in eV
286     aEckov [i] =1e-9*eV;    //Ckov energy in GeV
287     aAbsRad[i]=pRaAF->Eval(eV); aIdxRad[i]=1.292;//pRaIF->Eval(eV,20);      //Simulation for 20 degress C       
288     aAbsWin[i]=pWiAF->Eval(eV); aIdxWin[i]=1.5787;//pWiIF->Eval(eV);
289     aAbsGap[i]=pGaAF->Eval(eV); aIdxGap[i]=1.0005;//pGaIF->Eval(eV);    aQeAll[i] =1;                     //QE for all other materials except for PC must be 1.  
290     aAbsMet[i] =0.0001;                aIdxMet[i]=0;                                             //metal ref idx must be 0 in order to reflect photon
291                                        aIdxPc [i]=1;           aQePc [i]=pQeF->Eval(eV);         //PC ref idx must be 1 in order to apply photon to QE conversion 
292                                        
293   }
294   gMC->SetCerenkov((*fIdtmed)[kC6F14]    , kNbins, aEckov, aAbsRad  , aQeAll , aIdxRad );    
295   gMC->SetCerenkov((*fIdtmed)[kSiO2]     , kNbins, aEckov, aAbsWin  , aQeAll , aIdxWin );    
296   gMC->SetCerenkov((*fIdtmed)[kCH4]      , kNbins, aEckov, aAbsGap  , aQeAll , aIdxGap );    
297   gMC->SetCerenkov((*fIdtmed)[kCu]       , kNbins, aEckov, aAbsMet  , aQeAll , aIdxMet );    
298   gMC->SetCerenkov((*fIdtmed)[kW]        , kNbins, aEckov, aAbsMet  , aQeAll , aIdxMet ); //n=0 means reflect photons       
299   gMC->SetCerenkov((*fIdtmed)[kCsI]      , kNbins, aEckov, aAbsMet  , aQePc  , aIdxPc  ); //n=1 means convert photons    
300   gMC->SetCerenkov((*fIdtmed)[kAl]       , kNbins, aEckov, aAbsMet  , aQeAll , aIdxMet );    
301   delete pRaAF;delete pWiAF;delete pGaAF; delete pRaIF; delete pWiIF; delete pGaIF; delete pQeF;
302
303   AliDebug(1,"Stop v1 HMPID.");    
304 }
305 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
306 Bool_t AliHMPIDv1::IsLostByFresnel()
307 {
308 // Calculate probability for the photon to be lost by Fresnel reflection.
309   TLorentzVector p4;
310   Double_t mom[3],localMom[3];
311   gMC->TrackMomentum(p4);   mom[0]=p4(1);   mom[1]=p4(2);   mom[2]=p4(3);
312   localMom[0]=0; localMom[1]=0; localMom[2]=0;
313   gMC->Gmtod(mom,localMom,2);
314   Double_t localTc    = localMom[0]*localMom[0]+localMom[2]*localMom[2];
315   Double_t localTheta = TMath::ATan2(TMath::Sqrt(localTc),localMom[1]);
316   Double_t cotheta = TMath::Abs(TMath::Cos(localTheta));
317   if(gMC->GetRandom()->Rndm() < Fresnel(p4.E()*1e9,cotheta,1)){
318     AliDebug(1,"Photon lost");
319     return kTRUE;
320   }else
321     return kFALSE;
322 }//IsLostByFresnel()
323 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
324 void AliHMPIDv1::GenFee(Int_t iCh,Float_t eloss)
325 {
326 // Generate FeedBack photons for the current particle. To be invoked from StepManager().
327 // eloss=0 means photon so only pulse height distribution is to be analysed. This one is done in AliHMPIDParam::TotQdc()
328   TLorentzVector x4;
329   gMC->TrackPosition(x4); 
330   Int_t iNphotons=gMC->GetRandom()->Poisson(0.02*200); eloss++; iCh++;   //??????????????????????
331   AliDebug(1,Form("N photons=%i",iNphotons));
332   Int_t j;
333   Float_t cthf, phif, enfp = 0, sthf, e1[3], e2[3], e3[3], vmod, uswop,dir[3], phi,pol[3], mom[4];
334 //Generate photons
335   for(Int_t i=0;i<iNphotons;i++){//feedbacks loop
336     Double_t ranf[2];
337     gMC->GetRandom()->RndmArray(2,ranf);    //Sample direction
338     cthf=ranf[0]*2-1.0;
339     if(cthf<0) continue;
340     sthf = TMath::Sqrt((1 - cthf) * (1 + cthf));
341     phif = ranf[1] * 2 * TMath::Pi();
342     
343     if(Double_t randomNumber=gMC->GetRandom()->Rndm()<=0.57)
344       enfp = 7.5e-9;
345     else if(randomNumber<=0.7)
346       enfp = 6.4e-9;
347     else
348       enfp = 7.9e-9;
349     
350
351     dir[0] = sthf * TMath::Sin(phif);    dir[1] = cthf;    dir[2] = sthf * TMath::Cos(phif);
352     gMC->Gdtom(dir, mom, 2);
353     mom[0]*=enfp;    mom[1]*=enfp;    mom[2]*=enfp;
354     mom[3] = TMath::Sqrt(mom[0]*mom[0]+mom[1]*mom[1]+mom[2]*mom[2]);
355     
356     // Polarisation
357     e1[0]=      0;    e1[1]=-dir[2];    e1[2]= dir[1];
358     e2[0]=-dir[1];    e2[1]= dir[0];    e2[2]=      0;
359     e3[0]= dir[1];    e3[1]=      0;    e3[2]=-dir[0];
360     
361     vmod=0;
362     for(j=0;j<3;j++) vmod+=e1[j]*e1[j];
363     if (!vmod) for(j=0;j<3;j++) {
364       uswop=e1[j];
365       e1[j]=e3[j];
366       e3[j]=uswop;
367     }
368     vmod=0;
369     for(j=0;j<3;j++) vmod+=e2[j]*e2[j];
370     if (!vmod) for(j=0;j<3;j++) {
371       uswop=e2[j];
372       e2[j]=e3[j];
373       e3[j]=uswop;
374     }
375     
376     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;    
377     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;
378     
379     phi = gMC->GetRandom()->Rndm()* 2 * TMath::Pi();
380     for(j=0;j<3;j++) pol[j]=e1[j]*TMath::Sin(phi)+e2[j]*TMath::Cos(phi);
381     gMC->Gdtom(pol, pol, 2);
382     Int_t outputNtracksStored;    
383     gAlice->GetMCApp()->PushTrack(1,                             //transport
384                      gAlice->GetMCApp()->GetCurrentTrackNumber(),//parent track 
385                      kFeedback,                                  //PID
386                      mom[0],mom[1],mom[2],mom[3],                //track momentum  
387                      x4.X(),x4.Y(),x4.Z(),x4.T(),                //track origin 
388                      pol[0],pol[1],pol[2],                       //polarization
389                      kPFeedBackPhoton,                           //process ID   
390                      outputNtracksStored,                        //on return how many new photons stored on stack
391                      1.0);                                       //weight
392   }//feedbacks loop
393   AliDebug(1,"Stop.");
394 }//GenerateFeedbacks()
395 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
396 void AliHMPIDv1::Hits2SDigits()
397 {
398 // Interface methode ivoked from AliSimulation to create a list of sdigits corresponding to list of hits. Every hit generates one or more sdigits.
399 // Arguments: none
400 //   Returns: none   
401   AliDebug(1,"Start.");
402   for(Int_t iEvt=0;iEvt < GetLoader()->GetRunLoader()->GetNumberOfEvents();iEvt++){                //events loop
403     GetLoader()->GetRunLoader()->GetEvent(iEvt);                          //get next event
404   
405     if(!GetLoader()->TreeH()) {GetLoader()->LoadHits();                    }
406     if(!GetLoader()->TreeS()) {GetLoader()->MakeTree("S"); MakeBranch("S");}//to
407           
408     for(Int_t iPrimN=0;iPrimN<GetLoader()->TreeH()->GetEntries();iPrimN++){//prims loop
409       GetLoader()->TreeH()->GetEntry(iPrimN);
410       Hit2Sdi(Hits(),SdiLst());
411     }//prims loop
412     GetLoader()->TreeS()->Fill();
413     GetLoader()->WriteSDigits("OVERWRITE");
414     SdiReset();
415   }//events loop  
416   GetLoader()->UnloadHits();
417   GetLoader()->UnloadSDigits();  
418   AliDebug(1,"Stop.");
419 }//Hits2SDigits()
420 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
421 void AliHMPIDv1::Hit2Sdi(TClonesArray *pHitLst,TClonesArray *pSdiLst)
422 {
423 // Converts list of hits to list of sdigits.  For each hit in a loop the following steps are done:
424 // - calcultion of the total charge induced by the hit
425 // - determination of the pad contaning the hit and shifting hit y position to the nearest anod wire y
426 // - defining a set of pads affected (up to 9 including the hitted pad)    
427 // - calculating charge induced to all those pads using integrated Mathieson distribution and creating sdigit  
428 // Arguments: pHitLst  - list of hits provided not empty
429 //            pSDigLst - list of sdigits where to store the results
430 //   Returns: none         
431   for(Int_t iHit=0;iHit<pHitLst->GetEntries();iHit++){         //hits loop
432     AliHMPIDHit *pHit=(AliHMPIDHit*)pHitLst->At(iHit);           //get pointer to current hit   
433     AliHMPIDDigit::Hit2Sdi(pHit,pSdiLst);                       //convert this hit to list of sdigits     
434   }//hits loop loop
435 }//Hits2SDigs() for TVector2
436 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
437 void AliHMPIDv1::Digits2Raw()
438 {
439 // Creates raw data files in DDL format. Invoked by AliSimulation where loop over events is done
440 // Arguments: none
441 //   Returns: none    
442   AliDebug(1,"Start.");
443   GetLoader()->LoadDigits();
444   TTree * treeD = GetLoader()->TreeD();
445   if(!treeD) {
446     AliError("No digits tree!");
447     return;
448   }
449   treeD->GetEntry(0);
450   
451   ofstream file[AliHMPIDDigit::kNddls];   //output streams
452   Int_t    cnt[AliHMPIDDigit::kNddls];        //data words counters for DDLs
453   AliRawDataHeader header; //empty DDL header
454   UInt_t w32=0;            //32 bits data word 
455   
456   for(Int_t i=0;i<AliHMPIDDigit::kNddls;i++){        
457     file[i].open(AliDAQ::DdlFileName(GetName(),i)); //open all 14 DDL in parallel
458     file[i].write((char*)&header,sizeof(header));   //write dummy header as place holder, actual will be written later when total size of DDL is known
459     cnt[i]=0;                                       //reset counters
460   }
461   
462   for(Int_t iCh=0;iCh<7;iCh++)
463     for(Int_t iDig=0;iDig<DigLst(iCh)->GetEntriesFast();iDig++){//digits loop for a given chamber
464       AliHMPIDDigit *pDig=(AliHMPIDDigit*)DigLst(iCh)->At(iDig);
465       Int_t ddl=pDig->Raw(w32);                             //ddl is 0..13 
466       file[ddl].write((char*)&w32,sizeof(w32));  cnt[ddl]++;//write formated digit to the propriate file (as decided in Dig2Raw) and increment corresponding counter
467     }//digits 
468     
469     
470   for(Int_t i=0;i<AliHMPIDDigit::kNddls;i++){
471     header.fSize=sizeof(header)+cnt[i]*sizeof(w32);                //now calculate total number of bytes for each DDL file  
472     header.SetAttribute(0); 
473     file[i].seekp(0); file[i].write((char*)&header,sizeof(header));//rewrite DDL header with fSize field properly set
474     file[i].close();                                               //close DDL file
475   }
476   GetLoader()->UnloadDigits();
477   AliDebug(1,"Stop.");      
478 }//Digits2Raw()
479 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
480 Float_t AliHMPIDv1::Fresnel(Float_t ene,Float_t pdoti, Bool_t pola)
481 {
482 // Correction for Fresnel   ???????????
483 // Arguments:   ene - photon energy [GeV],
484 //              PDOTI=COS(INC.ANG.), PDOTR=COS(POL.PLANE ROT.ANG.)
485 //   Returns:  
486     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,
487                       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,
488                       7.8,7.9,8.0,8.1,8.2,8.3,8.4,8.5};
489     Float_t csin[36] = {2.14,2.21,2.33,2.48,2.76,2.97,2.99,2.59,2.81,3.05,
490                         2.86,2.53,2.55,2.66,2.79,2.96,3.18,3.05,2.84,2.81,2.38,2.11,
491                         2.01,2.13,2.39,2.73,3.08,3.15,2.95,2.73,2.56,2.41,2.12,1.95,
492                         1.72,1.53};
493     Float_t csik[36] = {0.,0.,0.,0.,0.,0.196,0.408,0.208,0.118,0.49,0.784,0.543,
494                         0.424,0.404,0.371,0.514,0.922,1.102,1.139,1.376,1.461,1.253,0.878,
495                         0.69,0.612,0.649,0.824,1.347,1.571,1.678,1.763,1.857,1.824,1.824,
496                         1.714,1.498};
497     Float_t xe=ene;
498     Int_t  j=Int_t(xe*10)-49;
499     Float_t cn=csin[j]+((csin[j+1]-csin[j])/0.1)*(xe-en[j]);
500     Float_t ck=csik[j]+((csik[j+1]-csik[j])/0.1)*(xe-en[j]);
501
502     //FORMULAE FROM HANDBOOK OF OPTICS, 33.23 OR
503     //W.R. HUNTER, J.O.S.A. 54 (1964),15 , J.O.S.A. 55(1965),1197
504
505     Float_t sinin=TMath::Sqrt(1-pdoti*pdoti);
506     Float_t tanin=sinin/pdoti;
507
508     Float_t c1=cn*cn-ck*ck-sinin*sinin;
509     Float_t c2=4*cn*cn*ck*ck;
510     Float_t aO=TMath::Sqrt(0.5*(TMath::Sqrt(c1*c1+c2)+c1));
511     Float_t b2=0.5*(TMath::Sqrt(c1*c1+c2)-c1);
512     
513     Float_t rs=((aO-pdoti)*(aO-pdoti)+b2)/((aO+pdoti)*(aO+pdoti)+b2);
514     Float_t rp=rs*((aO-sinin*tanin)*(aO-sinin*tanin)+b2)/((aO+sinin*tanin)*(aO+sinin*tanin)+b2);
515     
516
517     //CORRECTION FACTOR FOR SURFACE ROUGHNESS
518     //B.J. STAGG  APPLIED OPTICS, 30(1991),4113
519
520     Float_t sigraf=18.;
521     Float_t lamb=1240/ene;
522     Float_t fresn;
523  
524     Float_t  rO=TMath::Exp(-(4*TMath::Pi()*pdoti*sigraf/lamb)*(4*TMath::Pi()*pdoti*sigraf/lamb));
525
526     if(pola)
527     {
528         Float_t pdotr=0.8;                                 //DEGREE OF POLARIZATION : 1->P , -1->S
529         fresn=0.5*(rp*(1+pdotr)+rs*(1-pdotr));
530     }
531     else
532         fresn=0.5*(rp+rs);
533       
534     fresn = fresn*rO;
535     return fresn;
536 }//Fresnel()
537 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
538 void AliHMPIDv1::Print(Option_t *option)const
539 {
540 // Debug printout
541   TObject::Print(option);
542 }//void AliHMPID::Print(Option_t *option)const
543 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
544 Bool_t AliHMPIDv1::Raw2SDigits(AliRawReader *pRR)
545 {
546 // Interface methode ivoked from AliSimulation to create a list of sdigits from raw digits. Events loop is done in AliSimulation
547 // Arguments: pRR- raw reader 
548 //   Returns: kTRUE on success (currently ignored in AliSimulation::ConvertRaw2SDigits())      
549   AliHMPIDDigit sdi; //tmp sdigit, raw digit will be converted to it
550   
551   if(!GetLoader()->TreeS()) {MakeTree("S");  MakeBranch("S");}
552     
553   TClonesArray *pSdiLst=SdiLst(); Int_t iSdiCnt=0; //tmp list of sdigits for all chambers
554   pRR->Select("HMPID",0,13);//select all HMPID DDL files
555   UInt_t w32=0;
556   while(pRR->ReadNextInt(w32)){//raw records loop (in selected DDL files)
557     UInt_t ddl=pRR->GetDDLID(); //returns 0,1,2 ... 13
558     sdi.Raw(ddl,w32);  
559     new((*pSdiLst)[iSdiCnt++]) AliHMPIDDigit(sdi); //add this digit to the tmp list
560   }//raw records loop
561   GetLoader()->TreeS()->Fill(); GetLoader()->WriteSDigits("OVERWRITE");//write out sdigits
562   SdiReset();
563   return kTRUE;
564 }//Raw2SDigits
565 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
566 void AliHMPIDv1::StepCount()
567 {
568 // Count number of ckovs created  
569 }
570 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
571 void AliHMPIDv1::StepHistory()
572 {
573 // This methode is invoked from StepManager() in order to print out 
574   static Int_t iStepN;
575   const char *sParticle;
576   switch(gMC->TrackPid()){
577     case kProton:      sParticle="PROTON"    ;break;
578     case kNeutron:     sParticle="neutron"   ;break;
579     case kGamma:       sParticle="gamma"     ;break;
580     case kCerenkov:    sParticle="CKOV"    ;break;
581     case kPi0:         sParticle="Pi0"       ;break;  
582     case kPiPlus:      sParticle="Pi+"       ;break;  
583     case kPiMinus:     sParticle="Pi-"       ;break;  
584     case kElectron:    sParticle="electron"  ;break;  
585     default:           sParticle="not known" ;break;
586   }
587
588   TString flag="fanny combination";
589   if(gMC->IsTrackAlive())
590       if(gMC->IsTrackEntering())      flag="enters to";
591       else if(gMC->IsTrackExiting())  flag="exits from";
592       else if(gMC->IsTrackInside())   flag="inside";
593   else
594       if(gMC->IsTrackStop())          flag="stoped in";        
595   
596   Int_t vid=0,copy=0;
597   TString path=gMC->CurrentVolName(); path.Prepend("-");path.Prepend(gMC->CurrentVolOffName(1));//current volume and his mother are always there
598   vid=gMC->CurrentVolOffID(2,copy);  if(vid) {path.Prepend("-");path.Prepend(gMC->VolName(vid));}
599   vid=gMC->CurrentVolOffID(3,copy);  if(vid) {path.Prepend("-");path.Prepend(gMC->VolName(vid));}
600   
601   Printf("Step %i: %s (%i) %s %s m=%.6f GeV q=%.1f dEdX=%.4f",iStepN,sParticle,gMC->TrackPid(),flag.Data(),path.Data(),gMC->TrackMass(),gMC->TrackCharge(),gMC->Edep()*1e9);
602   
603   Printf("Step %i: tid=%i flags alive=%i disap=%i enter=%i exit=%i inside=%i out=%i stop=%i new=%i",
604                             iStepN, gAlice->GetMCApp()->GetCurrentTrackNumber(),
605                             gMC->IsTrackAlive(), gMC->IsTrackDisappeared(),gMC->IsTrackEntering(), gMC->IsTrackExiting(),
606                             gMC->IsTrackInside(),gMC->IsTrackOut(),        gMC->IsTrackStop(),     gMC->IsNewTrack());
607   
608   Float_t a,z,den,rad,abs; a=z=den=rad=abs=-1;
609   Int_t mid=gMC->CurrentMaterial(a,z,den,rad,abs);
610   Printf("Step %i: id=%i a=%7.2f z=%7.2f den=%9.4f rad=%9.2f abs=%9.2f\n\n",iStepN,mid,a,z,den,rad,abs);
611   iStepN++;
612 }//StepHistory()
613 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
614 void AliHMPIDv1::StepManager()
615 {
616 // Full Step Manager.
617 // Arguments: none
618 //   Returns: none           
619 //  StepHistory(); return; //uncomment to print tracks history
620 //  StepCount(); return;     //uncomment to count photons
621   
622   Int_t   copy; //volume copy aka node
623   
624 //Treat photons    
625   if((gMC->TrackPid()==kCerenkov||gMC->TrackPid()==kFeedback)&&gMC->CurrentVolID(copy)==fIdPc){  //photon (Ckov or feedback) hit PC (fIdPc)
626     if(gMC->Edep()>0){                                                                           //photon survided QE test i.e. produces electron
627       if(IsLostByFresnel()){ gMC->StopTrack(); return;}                                          //photon lost due to fersnel reflection on PC       
628                        gMC->CurrentVolOffID(2,copy);                                             //current chamber since geomtry tree is HMPID-Rppf-Rpc
629       Int_t   tid=     gMC->GetStack()->GetCurrentTrackNumber();                                 //take TID
630       Int_t   pid=     gMC->TrackPid();                                                          //take PID
631       Float_t etot=    gMC->Etot();                                                              //total hpoton energy, [GeV] 
632       Double_t x[3];   gMC->TrackPosition(x[0],x[1],x[2]);                                       //take MARS position at entrance to PC
633       Float_t xl,yl;   AliHMPIDParam::Instance()->Mars2Lors(copy,x,xl,yl);                        //take LORS position
634       new((*fHits)[fNhits++])AliHMPIDHit(copy,etot,pid,tid,xl,yl,x);                              //HIT for photon, position at PC
635       GenFee(copy);                                                                              //generate feedback photons
636     }//photon hit PC and DE >0 
637   }//photon hit PC
638   
639 //Treat charged particles  
640   static Double_t dEdX;                                                                           //need to store mip parameters between different steps    
641   static Double_t in[3];
642   if(gMC->TrackCharge() && gMC->CurrentVolID(copy)==fIdAmpGap){                                   //charged particle in amplification gap (fIdAmpGap)
643     if(gMC->IsTrackEntering()||gMC->IsNewTrack()) {                                               //entering or newly created
644       dEdX=0;                                                                                     //reset dEdX collector                         
645       gMC->TrackPosition(in[0],in[1],in[2]);                                                      //take position at the entrance
646     }else if(gMC->IsTrackExiting()||gMC->IsTrackStop()||gMC->IsTrackDisappeared()){               //exiting or disappeared
647       dEdX              +=gMC->Edep();                                                            //take into account last step dEdX
648                           gMC->CurrentVolOffID(1,copy);                                           //take current chamber since geometry tree is HMPID-Rgap
649       Int_t tid=          gMC->GetStack()->GetCurrentTrackNumber();                               //take TID
650       Int_t pid=          gMC->TrackPid();                                                        //take PID
651       Double_t out[3];    gMC->TrackPosition(out[0],out[1],out[2]);                               //take MARS position at exit
652       out[0]=0.5*(out[0]+in[0]); out[1]=0.5*(out[1]+in[1]); out[1]=0.5*(out[1]+in[1]);            //take hit position at the anod plane
653       Float_t xl,yl;AliHMPIDParam::Instance()->Mars2Lors(copy,out,xl,yl);                          //take LORS position
654       new((*fHits)[fNhits++])AliHMPIDHit(copy,dEdX,pid,tid,xl,yl,out);                             //HIT for MIP, position near anod plane 
655       GenFee(copy,dEdX);                                                                          //generate feedback photons
656     }else                                                                                         //just going inside
657       dEdX          += gMC->Edep();                                                               //collect this step dEdX 
658   }//MIP in GAP
659 }//StepManager()
660 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++