]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDv3.cxx
Removing centrality check from default:
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDv3.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 "AliHMPIDv3.h"         //class header
18 #include "AliHMPIDParam.h"      //StepManager()
19 #include "AliHMPIDHit.h"        //Hits2SDigs(),StepManager()
20 #include "AliHMPIDDigit.h"      //Digits2Raw(), Raw2SDigits()
21 #include "AliHMPIDRawStream.h"  //Digits2Raw(), Raw2SDigits()
22 #include "AliRawReader.h"       //Raw2SDigits()
23 #include "AliTrackReference.h"
24 #include <TVirtualMC.h>         //StepManager() for TVirtualMC::GetMC()
25 #include <TPDGCode.h>           //StepHistory() 
26 #include <AliStack.h>           //StepManager(),Hits2SDigits()78.6
27 #include <AliLoader.h>          //Hits2SDigits()
28 #include <AliRunLoader.h>       //Hits2SDigits()
29 #include <AliMC.h>              //StepManager()      
30 #include <AliRun.h>             //CreateMaterials()    
31 #include <AliMagF.h>            //CreateMaterials()
32 #include "AliGeomManager.h"   //AddAlignableVolumes()
33 #include <AliCDBEntry.h>        //CreateMaterials()
34 #include <AliCDBManager.h>      //CreateMaterials()
35 #include <TF1.h>                //DefineOpticalProperties()
36 #include <TF2.h>                //DefineOpticalProperties()
37 #include <TGeoCompositeShape.h> //CradleBaseVolume()
38 #include <TGeoGlobalMagField.h>
39 #include <TGeoPhysicalNode.h>   //AddAlignableVolumes()
40 #include <TGeoXtru.h>           //CradleBaseVolume()
41 #include <TLorentzVector.h>     //IsLostByFresnel() 
42 #include <TString.h>            //StepManager()
43 #include <TTree.h>
44
45 ClassImp(AliHMPIDv3)    
46 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
47 void AliHMPIDv3::AddAlignableVolumes()const
48 {
49 // Associates the symbolic volume name with the corresponding volume path. Interface method from AliModule invoked from AliMC
50 // Arguments: none
51 //   Returns: none   
52   
53   AliGeomManager::ELayerID idHMPID = AliGeomManager::kHMPID;
54   Int_t modUID, modnum = 0;
55
56   TGeoHMatrix *pGm = new TGeoHMatrix;
57   Double_t trans[3]={0.5*131.24,0.5*126.16,0};                            //translation from LORS to TGeo RS (half size AllX, half size allY,0)
58   pGm->SetTranslation(trans);
59  
60   Double_t ph[7]={10.,10., 30.,30.,30. ,50.,50};
61
62   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++) {
63     modUID = AliGeomManager::LayerToVolUID(idHMPID,modnum++);
64     if(!gGeoManager->SetAlignableEntry(Form("/HMPID/Chamber%i",iCh),Form("ALIC_1/Hmp%i_0",iCh),modUID))
65             AliError("AliHMPIDv3::Unable to set alignable entry!!");  //aligment without AliCluster3D
66     //Get Tracking To Local matricies for alignment with AliCluster3D
67     TGeoPNEntry *eCh = gGeoManager->GetAlignableEntryByUID(modUID);
68     TGeoHMatrix *globMatrix = eCh->GetGlobalOrig();
69
70     //Double_t phi = 20.0 * ((iCh+1) / 3) + 10.0;
71     Double_t phi = ph[iCh];
72     TGeoHMatrix *t2l  = new TGeoHMatrix();
73     t2l->RotateZ(phi);
74     t2l->MultiplyLeft(&(globMatrix->Inverse()));
75     eCh->SetMatrix(t2l);
76   }//iCh loop
77   
78 }
79 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
80 void AliHMPIDv3::CreateMaterials()
81 {
82 // Definition of available HMPID materials  
83 // Arguments: none
84 //   Returns: none    
85   AliDebug(1,"Start v2 HMPID.");
86     
87     //clm update material definition later on from Antonello
88     
89 //data from PDG booklet 2002     density [gr/cm^3] rad len [cm] abs len [cm]    
90   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
91   Float_t aC6F14[2]={ 12.01 , 18.99} , zC6F14[2]={ 6 , 9}   , wC6F14[2]={6 , 14} , dC6F14=1.68    ; Int_t nC6F14=-2;
92   Float_t  aSiO2[2]={ 28.09 , 15.99} ,  zSiO2[2]={14 , 8}   ,  wSiO2[2]={1 ,  2} ,  dSiO2=2.64    ; Int_t  nSiO2=-2; 
93   Float_t   aCH4[2]={ 12.01 ,  1.01} ,   zCH4[2]={ 6 , 1}   ,   wCH4[2]={1 ,  4} ,   dCH4=7.17e-4 ; Int_t   nCH4=-2; 
94 // not necessary...PCB properties instead! Float_t   aCsI[2]={132.90 ,126.90} ,   zCsI[2]={55 ,53}   ,   wCsI[2]={1 ,  1} ,   dCsI=0.1     ; Int_t   nCsI=-2; 
95   
96   Float_t     aRoha = 12.01 ,   zRoha =  6 ,  dRoha =  0.10    ,   radRoha = 18.80 , absRoha =  86.3/dRoha; //special material- quasi quartz
97   Float_t       aCu = 63.55 ,   zCu   = 29 ,  dCu   =  8.96    ,   radCu   =  1.43 , absCu   = 134.9/dCu  ;
98   Float_t        aW =183.84 ,   zW    = 74 ,  dW    = 19.30    ,   radW    =  0.35 , absW    = 185.0/dW   ;
99   Float_t       aAl = 26.98 ,   zAl   = 13 ,  dAl   =  2.70    ,   radAl   =  8.90 , absAl   = 106.4/dAl  ;
100   Float_t       aAr = 39.94 ,   zAr   = 18 ,  dAr   =  1.396e-3,   radAr   =  14.0 , absAr   = 117.2/dAr  ;   
101
102     Int_t   matId=0;                           //tmp material id number
103     Int_t   unsens =  0, sens=1;               //sensitive or unsensitive medium
104     Int_t   itgfld = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Integ(); //type of field intergration 0 no field -1 user in guswim 1 Runge Kutta 2 helix 3 const field along z
105     Float_t maxfld = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Max();   //max field value
106     Float_t tmaxfd = -10.0;                    //max deflection angle due to magnetic field in one step
107     Float_t deemax = - 0.2;                    //max fractional energy loss in one step   
108     Float_t stemax = - 0.1;                    //max step allowed [cm]
109     Float_t epsil  =   0.001;                  //abs tracking precision [cm]   
110     Float_t stmin  = - 0.001;                  //min step size [cm] in continius process transport, negative value: choose it automatically
111
112     // PCB copmposed mainly by G10 (Si,C,H,O) -> CsI is negligible (<500nm thick)
113     // So what is called CsI has the optical properties of CsI, but the composition of G-10 (for delta elec, etc production...)
114     
115     Float_t aG10[4] = {28.09,12.01,1.01,16.00};
116     Float_t zG10[4] = {14.,  6.,  1.,  8.};
117     Float_t wG10[4] = {0.129060,0.515016,0.061873,0.294050};
118     Float_t dG10    = 1.7;
119     Int_t   nG10    = 4;
120     
121     AliMixture(++matId,"Air"  ,aAir  ,zAir  ,dAir  ,nAir  ,wAir  ); AliMedium(kAir  ,"Air"  ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
122     AliMixture(++matId,"C6F14",aC6F14,zC6F14,dC6F14,nC6F14,wC6F14); AliMedium(kC6F14,"C6F14",matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);      
123     AliMixture(++matId,"SiO2" ,aSiO2 ,zSiO2 ,dSiO2 ,nSiO2 ,wSiO2 ); AliMedium(kSiO2 ,"SiO2" ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);    
124     AliMixture(++matId,"CH4"  ,aCH4  ,zCH4  ,dCH4  ,nCH4  ,wCH4  ); AliMedium(kCH4  ,"CH4"  ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);  
125 //    AliMixture(++matId,"CsI"  ,aCsI  ,zCsI  ,dCsI  ,nCsI  ,wCsI  ); AliMedium(kCsI  ,"CsI"  ,matId,   sens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);//sensitive
126     AliMixture(++matId,"CsI+PCB",aG10  , zG10, dG10,nG10   ,wG10   ); AliMedium(kCsI  ,"CsI"  ,matId,   sens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);//sensitive
127
128     AliMixture(++matId ,"Neo" ,aSiO2 ,zSiO2 ,dSiO2 ,nSiO2 ,wSiO2 ); AliMedium(kNeo  ,"Neo"  ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin); //clm neoceram
129     AliMaterial(++matId,"Roha",aRoha,zRoha,dRoha,radRoha,absRoha);  AliMedium(kRoha ,"Roha" ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin); //Roha->honeycomb
130
131
132     AliMaterial(++matId,"Cu"  ,aCu  ,zCu  ,dCu  ,radCu  ,absCu  );  AliMedium(kCu  ,"Cu"  , matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
133     AliMaterial(++matId,"W"   ,aW   ,zW   ,dW   ,radW   ,absW   );  AliMedium(kW   ,"W"   , matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
134     AliMaterial(++matId,"Al"  ,aAl  ,zAl  ,dAl  ,radAl  ,absAl  );  AliMedium(kAl  ,"Al"  , matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
135     AliMaterial(++matId,"Ar"  ,aAr  ,zAr  ,dAr  ,radAr  ,absAr  );  AliMedium(kAr  ,"Ar"  , matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
136         
137 }//void AliHMPID::CreateMaterials()
138 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
139 //void AliHMPIDv3::InitProperties()
140 //{
141 /*
142 * HMPID
143 * ====
144 *
145 *       GAM   ELEC  NHAD   CHAD  MUON  EBREM MUHAB  EDEL  MUDEL MUPA ANNI BREM COMP DCAY DRAY HADR LOSS MULS PAIR PHOT RAYL
146 * Quarz Window        (>1000 keV delta-electrons)
147 HMPID  3  1.e-4 1.e-4 1.e-4  -1.   1.e-4 -1.   -1.    1.e-3 1.e-3 -1.  -1   -1   -1   -1   1    -1   1    -1   -1   -1   -1 
148 * Freon Radiator      (>  500 keV delta-electrons)
149 HMPID  4  1.e-4 1.e-4 1.e-4  -1.   1.e-4 -1.   -1.    5.e-4 5.e-4 -1.  -1   -1   -1   -1   1    -1   1    -1   -1   -1   -1 
150 * Methane Gap         (>  100 keV delta-electrons)
151 HMPID  5  5.e-5 1.e-5 1.e-4 -1.   1.e-4 -1.   -1.     1.e-4 1.e-4 -1.  -1   -1   -1   -1   1    -1   1    -1   -1   -1   -1 
152 * Sensitive Volume    (>  50 keV delta-electrons)
153 HMPID  9  1.e-5 1.e-5 1.e-4  -1.   1.e-4 -1.   -1.    5.e-5 5.e-5 -1.  -1   -1   -1   -1   1    -1   1    -1   -1   -1   -1 
154 * CSI    (>  50 keV delta-electrons)
155 HMPID  6  1.e-5 1.e-5 1.e-4  -1.   1.e-4 -1.   -1.    5.e-5 5.e-5 -1.  -1   -1   -1   -1   1    -1   1    -1   -1   -1   -1 
156 * PCB backplane   (>  50 keV delta-electrons)
157 HMPID 12  1.e-5 1.e-5 1.e-4  -1.   1.e-4 -1.   -1.    5.e-5 5.e-5 -1.  -1   -1   -1   -1   1    -1   1    -1   -1   -1   -1 
158
159     Int_t *idtmed = fIdtmed->GetArray();
160     Int_t imed;
161     
162     imed = kSiO2;   // * Quarz Window        (>1000 keV delta-electrons)
163     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTGAM",1.e-4);
164     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTELE",1.e-4);
165     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTNEU",1.e-4);
166     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTMUO",1.e-4);    
167     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "DCUTE" ,1.e-3);    
168     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTHAD",1.e-3);    
169     
170     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "DRAY",1);    
171     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "LOSS",1);    
172
173     imed = kC6F14;  // * Freon Radiator      (>  500 keV delta-electrons)
174     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTGAM",1.e-4);
175     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTELE",1.e-4);
176     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTNEU",1.e-4);
177     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTMUO",1.e-4);    
178     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "DCUTE" ,5.e-4);    
179     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTHAD",5.e-4);    
180     
181     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "DRAY",1);    
182     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "LOSS",1);    
183     
184     imed = kCH4;  // * Methane Gap         (>  100 keV delta-electrons)
185     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTGAM",5.e-5);
186     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTELE",5.e-5);
187     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTNEU",1.e-4);
188     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTMUO",1.e-4);    
189     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "DCUTE" ,1.e-4);    
190     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTHAD",1.e-4);    
191     
192     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "DRAY",1);    
193     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "LOSS",1);    
194     
195     imed = kCsI;  // * CSI    (>  50 keV delta-electrons)
196     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTGAM",1.e-5);
197     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTELE",1.e-5);
198     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTNEU",1.e-4);
199     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTMUO",1.e-4);    
200     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "DCUTE" ,5.e-5);    
201     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTHAD",5.e-5);    
202     
203     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "DRAY",1);    
204     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "LOSS",1);
205     
206     imed = kAl;  // * Alluminium    (>  50 keV delta-electrons)
207     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTGAM",1.e-5);
208     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTELE",1.e-5);
209     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTNEU",1.e-4);
210     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTMUO",1.e-4);    
211     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "DCUTE" ,5.e-5);    
212     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTHAD",5.e-5);    
213     
214     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "DRAY",1);    
215     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "LOSS",1);    
216     
217     imed = kCu;  // * Copper       (>  50 keV delta-electrons)
218     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTGAM",1.e-5);
219     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTELE",1.e-5);
220     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTNEU",1.e-4);
221     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTMUO",1.e-4);    
222     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "DCUTE" ,5.e-5);    
223     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTHAD",5.e-5);    
224     
225     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "DRAY",1);    
226     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "LOSS",1);    
227     
228     imed = kW;  // * Tungsten     (>  50 keV delta-electrons)
229     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTGAM",1.e-5);
230     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTELE",1.e-5);
231     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTNEU",1.e-4);
232     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTMUO",1.e-4);    
233     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "DCUTE" ,5.e-5);    
234     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "CUTHAD",5.e-5);    
235     
236     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "DRAY",1);    
237     TVirtualMC::GetMC()->Gstpar(idtmed[imed], "LOSS",1);    
238     
239 }*/
240 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
241 void AliHMPIDv3::CreateGeometry()
242 {
243 //Creates detailed geometry simulation (currently GEANT volumes tree)         
244 //includind the HMPID cradle
245
246   AliDebug(1,"Start main.");
247   if(!TVirtualMC::GetMC()->IsRootGeometrySupported()) return;                
248
249   TGeoVolume *hmpcradle = CreateCradle();
250   TString title=GetTitle();
251   if(title.Contains("TestBeam")){
252    TGeoVolume *hmpid = CreateChamber(3);
253     gGeoManager->GetVolume("ALIC")->AddNode(hmpid,0);
254   }else{
255     for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){//place 7 chambers
256       TGeoVolume *hmpid = CreateChamber(iCh);
257       TGeoHMatrix *pMatrix=new TGeoHMatrix;
258       IdealPosition(iCh,pMatrix);
259       gGeoManager->GetVolume("ALIC")->AddNode(hmpid,0,pMatrix);
260       if(iCh==1 || iCh == 3 || iCh == 5){ 
261      TGeoHMatrix *pCradleMatrix=new TGeoHMatrix;
262       IdealPositionCradle(iCh,pCradleMatrix);
263       gGeoManager->GetVolume("ALIC")->AddNode(hmpcradle,iCh,pCradleMatrix);
264       }  
265     }
266   }
267   AliDebug(1,"Stop v3. HMPID option");
268
269   
270
271 TGeoVolume * AliHMPIDv3::CreateChamber(Int_t number)
272 {
273   //Single module geometry building
274    
275  Double_t cm=1,mm=0.1*cm,um=0.001*mm;//default is cm
276  
277   TGeoVolume *hmp = new TGeoVolumeAssembly(Form("Hmp%i",number));
278  
279   TGeoMedium *al   =gGeoManager->GetMedium("HMPID_Al");    
280   TGeoMedium *ch4  =gGeoManager->GetMedium("HMPID_CH4");    
281   TGeoMedium *roha =gGeoManager->GetMedium("HMPID_Roha");   
282   TGeoMedium *neoc =gGeoManager->GetMedium("HMPID_Neo");
283   TGeoMedium *c6f14=gGeoManager->GetMedium("HMPID_C6F14");  
284   TGeoMedium *sio2 =gGeoManager->GetMedium("HMPID_SiO2");   
285   TGeoMedium *cu   =gGeoManager->GetMedium("HMPID_Cu");     
286   TGeoMedium *w    =gGeoManager->GetMedium("HMPID_W");      
287   TGeoMedium *csi  =gGeoManager->GetMedium("HMPID_CsI");    
288   TGeoMedium *ar   =gGeoManager->GetMedium("HMPID_Ar");     
289   
290
291   TGeoRotation *rot=new TGeoRotation("HwireRot"); rot->RotateY(90); //rotate wires around Y to be along X (initially along Z)
292   TGeoVolume *sbo=gGeoManager->MakeBox ("Hsbo",ch4  , 1419*mm/2 , 1378.00*mm/2 ,   50.5*mm/2);//2072P1
293   TGeoVolume *cov=gGeoManager->MakeBox ("Hcov",al   , 1419*mm/2 , 1378.00*mm/2 ,    0.5*mm/2);  
294   TGeoVolume *hon=gGeoManager->MakeBox ("Hhon",roha , 1359*mm/2 , 1318.00*mm/2 ,   49.5*mm/2);  
295   TGeoVolume *rad=gGeoManager->MakeBox ("Hrad",c6f14, 1330*mm/2 ,  413.00*mm/2 ,   24.0*mm/2); //2011P1
296   TGeoVolume *neo=gGeoManager->MakeBox ("Hneo",neoc , 1330*mm/2 ,  413.00*mm/2 ,    4.0*mm/2); 
297   TGeoVolume *win=gGeoManager->MakeBox ("Hwin",sio2 , 1330*mm/2 ,  413.00*mm/2 ,    5.0*mm/2); 
298   TGeoVolume *si1=gGeoManager->MakeBox ("Hsi1",sio2 , 1330*mm/2 ,    5.00*mm/2 ,   15.0*mm/2);    
299   TGeoVolume *si2=gGeoManager->MakeBox ("Hsi2",neoc ,   10*mm/2 ,  403.00*mm/2 ,   15.0*mm/2);    
300   TGeoVolume *spa=gGeoManager->MakeTube("Hspa",sio2 ,    0*mm   ,    5.00*mm   ,   15.0*mm/2);         
301   TGeoVolume *fr4=gGeoManager->MakeBox ("Hfr4",ch4  , 1407*mm/2 , 1366.00*mm/2 ,   15.0*mm/2);//2043P1 
302   TGeoVolume *f4a=gGeoManager->MakeBox ("Hf4a",al   , 1407*mm/2 , 1366.00*mm/2 ,   10.0*mm/2); 
303   TGeoVolume *f4i=gGeoManager->MakeBox ("Hf4i",ch4  , 1323*mm/2 , 1296.00*mm/2 ,   10.0*mm/2); 
304   TGeoVolume *col=gGeoManager->MakeTube("Hcol",cu   ,    0*mm   ,  100.00*um   , 1323.0*mm/2);
305   TGeoVolume *sec=gGeoManager->MakeBox ("Hsec",ch4  ,  648*mm/2 ,  411.00*mm/2 ,   6.2*mm/2);//sec=gap 2099P1 (6.2 = 4.45 + 0.05 (1/2 diameter wire)+1.7)
306  
307   Double_t cellx=8.04*mm,celly=8.4*mm;  Int_t nPadX=80, nPadY=48; 
308   TGeoVolume *gap=gGeoManager->MakeBox ("Hgap",ch4  ,  cellx*nPadX/2 ,  celly*nPadY/2 ,    6.2*mm/2); //x=8.04*80 y=8.4*48 z=pad+pad-ano+marign 2006p1  
309   TGeoVolume *row=        gap->Divide  ("Hrow",2,nPadY,0,0);//along Y->48 rows
310   TGeoVolume *cel=        row->Divide  (Form("Hcel%i",number),1,nPadX,0,0);//along X->80 cells
311   TGeoVolume *cat=gGeoManager->MakeTube("Hcat",cu   ,    0.00*mm   ,   50.00*um   ,    cellx/2); 
312   TGeoVolume *ano=gGeoManager->MakeTube("Hano",w    ,    0.00*mm   ,   20.00*um   ,    cellx/2); 
313   TGeoVolume *pad=gGeoManager->MakeBox (Form("Hpad%i",number),csi  ,    7.54*mm/2 ,    7.90*mm/2 ,    1.7*mm/2); //2006P1 PCB material...     
314   TGeoVolume *fr1=gGeoManager->MakeBox ("Hfr1",al   , 1463*mm/2 , 1422.00*mm/2 ,   58.3*mm/2);//2040P1 and pad plane is excluded (62 - 2 - 17)
315   TGeoVolume *fr1up=gGeoManager->MakeBox ("Hfr1up",ch4,(1426.00-37.00)*mm/2 , (1385.00-37.00)*mm/2 ,    20.0*mm/2);//2040P1
316   
317   TGeoVolume *fr1upcard=gGeoManager->MakeBox ("Hfr1upcard",ch4,662.*mm/2., 425.*mm/2. ,19.0*mm/2);//needed to set the gassiplex
318
319   TGeoVolume *fr1perUpBig=gGeoManager->MakeBox ("Hfr1perUpBig",ch4,1389*mm/2,35*mm/2,10*mm/2);    
320   TGeoVolume *fr1perUpSma=gGeoManager->MakeBox ("Hfr1perUpSma",ch4,35*mm/2,(1385-37-2*35)*mm/2,10*mm/2);
321         TGeoVolume *fr1perDowBig=gGeoManager->MakeBox ("Hfr1perDowBig",ch4,1389*mm/2,46*mm/2,2.3*mm/2);    
322   TGeoVolume *fr1perDowSma=gGeoManager->MakeBox ("Hfr1perDowSma",ch4,46*mm/2,(1385-37-2*46)*mm/2,2.3*mm/2);
323         
324         TGeoVolume *ppf=gGeoManager->MakeBox ("Hppf",al   ,  648*mm/2 ,  411.00*mm/2 ,   38.3*mm/2);//2001P2
325   TGeoVolume *lar=gGeoManager->MakeBox ("Hlar",ar   ,  181*mm/2 ,   89.25*mm/2 ,   38.3*mm/2);//2001P2
326   TGeoVolume *smo=gGeoManager->MakeBox ("Hsmo",ar   ,  114*mm/2 ,   89.25*mm/2 ,   38.3*mm/2);//2001P2
327                 
328   TGeoVolume *cufoil = gGeoManager->MakeBox("Hcufoil", csi, 662.*mm/2., 425.*mm/2., 1.*mm/2.);//PCB foil at the back of the ppf with holes for GASSIPLEX
329   TGeoVolume *rect = gGeoManager->MakeBox("Hrect",ch4, 48*mm/2, 19*mm/2., 1*mm/2.);
330
331                 
332         TGeoVolume *fr3=   gGeoManager->MakeBox("Hfr3",          al,  1463*mm/2,  1422*mm/2,  34*mm/2);//2041P1
333    TGeoVolume *fr3up=    gGeoManager->MakeBox("Hfr3up",     ch4, 1323*mm/2,  1282*mm/2,  20*mm/2);//2041P1
334    TGeoVolume *fr3down=gGeoManager->MakeBox("Hfr3down", ch4, 1437*mm/2,  1370*mm/2,  14*mm/2);//2041P1
335
336
337   TGeoVolume *proxgap1 = gGeoManager->MakeBox("Hproxgap1",ch4,1407*mm/2 , 1366.00*mm/2 ,(9.-7.5)*mm/2.);//methane volume between quartz and fr4
338   TGeoVolume *proxgap2 = gGeoManager->MakeBox("Hproxgap2",ch4,1407*mm/2 , 1366.00*mm/2 ,(81.7-6.2-34.-9.-7.5)*mm/2.);//methane volume between fr4 and Hgap(tot height(81.7) - Hsec (6.2)   - proxygap2 (34) - upper bound of fr4 (9+7.5))
339   
340
341 // ^ Y   z=         z=-12mm      z=98.25mm               ALIC->7xHmp (virtual)-->1xHsbo (virtual) --->2xHcov (real) 2072P1
342 // |  ____________________________________                                    |                   |-->1xHhon (real) 2072P1
343 // | |   ______     ____          ______  |                                   |
344 //   |  |      |   |    |   *    |      | |                                   |->3xHrad (virtual) --->1xHneo (real) 2011P1
345 //   |  |50.5mm|   |24mm|   *    |45.5mm| |                                   |                   |-->1xHwin (real) 2011P1
346 //   |  |      |   |    |   *    |      | |                                   |                   |-->2xHsi1 (real) 2011P1
347 //   |  |      |   |____|   *    |______| |                                   |                   |-->2xHsi2 (real) 2011P1
348 //   |  |      |    ____    *     ______  |                                   |                   |->30xHspa (real) 2011P1
349 //   |  |      |   |    |   *    |      | |                                   |
350 //   |  |      |   |    |   *    |      | |                                   |->1xHfr4 (vitual) --->1xHf4a (real)---->1xHf4i(virtual) 2043P1 
351 //   |  |  sb  |   | rad|   *    |      | |                                   |                  |-->322xHcol (real) 2043P1
352 //   |  |      |   |____|   *    |______| |                                   |
353 //   |  |      |    ____    *     ______  |                                   |->1xHfr1 (real) --> 6xHppf(real) ---->8xHlar (virtual) 2001P1
354 //   |  |      |   |    |   *    |      | |                                   |                                     |--->8xHsmo (virtual) 2001P1     
355 //   |  |      |   |    |   *    |      | |                                   |               
356 //   |  |      |   |    |   *    |      | |                                   |-> 6xHgap (virtual) --->48xHrow (virtual) -->80xHcel (virtual) -->4xHcat (real) from p84 TDR 
357 //   |  |______|   |____|   *    |______| |                                                                                                  |-->2xHano (real) from p84 TDR                                  
358 //   |____________________________________|                                                                                                  |-->1xHpad (real) from p84 TDR 
359 //                                                       --->Z 
360   hmp->AddNode(sbo ,1,new TGeoTranslation(   0*mm,   0*mm, -73.75*mm));                     //p.84 TDR
361      sbo->AddNode(hon ,1,new TGeoTranslation(  0*mm,0*mm,      0*mm)); //2072P1
362      sbo->AddNode(cov ,1,new TGeoTranslation(  0*mm,0*mm,    +25*mm)); 
363      sbo->AddNode(cov ,2,new TGeoTranslation(  0*mm,0*mm,    -25*mm)); 
364   hmp->AddNode(rad,2,new TGeoTranslation(   0*mm,+434*mm, -12.00*mm)); 
365   hmp->AddNode(rad,1,new TGeoTranslation(   0*mm,   0*mm, -12.00*mm)); 
366   hmp->AddNode(rad,0,new TGeoTranslation(   0*mm,-434*mm, -12.00*mm)); 
367     rad->AddNode(neo,1,new TGeoTranslation(   0*mm,   0*mm, -10.0*mm));
368     rad->AddNode(win,1,new TGeoTranslation(   0*mm,   0*mm,   9.5*mm));
369     rad->AddNode(si1,1,new TGeoTranslation(   0*mm,-204*mm,  -0.5*mm)); rad->AddNode(si1,2,new TGeoTranslation(   0*mm,+204*mm,  -0.5*mm));
370     rad->AddNode(si2,1,new TGeoTranslation(-660*mm,   0*mm,  -0.5*mm)); rad->AddNode(si2,2,new TGeoTranslation(+660*mm,   0*mm,  -0.5*mm));
371     for(Int_t i=0;i<3;i++) for(Int_t j=0;j<10;j++) rad->AddNode(spa,10*i+j,new TGeoTranslation(-1330*mm/2+116*mm+j*122*mm,(i-1)*105*mm,-0.5*mm));
372   hmp->AddNode(fr4,1,new TGeoTranslation(   0*mm,   0*mm,   9.00*mm));                     //p.84 TDR
373   for(int i=1;i<=322;i++)  fr4->AddNode(col,i,new TGeoCombiTrans( 0*mm, -1296/2*mm+i*4*mm,-5*mm,rot)); //F4 2043P1
374                            fr4->AddNode(f4a,1,new TGeoTranslation(   0*mm,0*mm, 2.5*mm));    
375                                         f4a->AddNode(f4i,1,new TGeoTranslation(   0*mm,0*mm,   0*mm));
376   hmp->AddNode(sec,4,new TGeoTranslation(-335*mm,+433*mm,  78.6*mm)); hmp->AddNode(sec,5,new TGeoTranslation(+335*mm,+433*mm,  78.6*mm));
377   hmp->AddNode(sec,2,new TGeoTranslation(-335*mm,   0*mm,  78.6*mm)); hmp->AddNode(sec,3,new TGeoTranslation(+335*mm,   0*mm,  78.6*mm));
378   hmp->AddNode(sec,0,new TGeoTranslation(-335*mm,-433*mm,  78.6*mm)); hmp->AddNode(sec,1,new TGeoTranslation(+335*mm,-433*mm,  78.6*mm));
379     sec->AddNode(gap,1,new TGeoTranslation(0,0,0.*mm));
380       cel->AddNode(cat,1,new TGeoCombiTrans (0,  3.15*mm , -2.70*mm , rot)); //4 cathode wires
381       cel->AddNode(ano,1,new TGeoCombiTrans (0,  2.00*mm , -0.29*mm , rot)); //2 anod wires
382       cel->AddNode(cat,2,new TGeoCombiTrans (0,  1.05*mm , -2.70*mm , rot)); 
383       cel->AddNode(cat,3,new TGeoCombiTrans (0, -1.05*mm , -2.70*mm , rot)); 
384       cel->AddNode(ano,2,new TGeoCombiTrans (0, -2.00*mm , -0.29*mm , rot)); 
385       cel->AddNode(cat,4,new TGeoCombiTrans (0, -3.15*mm , -2.70*mm , rot));   
386       cel->AddNode(pad,1,new TGeoTranslation(0,  0.00*mm ,  2.25*mm));       //1 pad  
387             
388   hmp->AddNode(fr1,1,new TGeoTranslation(0.,0.,(80.+1.7)*mm+58.3*mm/2.));
389                 fr1->AddNode(fr1up,1,new TGeoTranslation(0.,0.,(58.3*mm-20.00*mm)/2.));
390                 
391                 fr1->AddNode(fr1perUpBig,0,new TGeoTranslation(0.,(1385-37-35)*mm/2.,(58.3*mm-20.00*2*mm-10.0*mm)/2.));
392                 fr1->AddNode(fr1perUpSma,0,new TGeoTranslation((1426-37-35)*mm/2.,0.,(58.3*mm-20.00*2*mm-10.0*mm)/2.));
393                 fr1->AddNode(fr1perUpBig,1,new TGeoTranslation(0.,-(1385-37-35)*mm/2.,(58.3*mm-20.00*2*mm-10.0*mm)/2.));
394                 fr1->AddNode(fr1perUpSma,1,new TGeoTranslation(-(1426-37-35)*mm/2.,0.,(58.3*mm-20.00*2*mm-10.0*mm)/2.));
395                 
396           fr1->AddNode(fr1perDowBig,0,new TGeoTranslation(0.,(1385-37)*mm/2.,(-58.3*mm+2.3*mm)/2.));
397                 fr1->AddNode(fr1perDowSma,0,new TGeoTranslation((1426-37)*mm/2.,0.,(-58.3*mm+2.3*mm)/2.));
398           fr1->AddNode(fr1perDowBig,1,new TGeoTranslation(0.,-(1385-37)*mm/2.,(-58.3*mm+2.3*mm)/2.));
399                 fr1->AddNode(fr1perDowSma,1,new TGeoTranslation(-(1426-37)*mm/2.,0.,(-58.3*mm+2.3*mm)/2.));             
400                         
401           fr1->AddNode(ppf,4,new TGeoTranslation(-335*mm,433*mm,(-58.3+38.3)*mm/2.));  fr1->AddNode(ppf,5,new TGeoTranslation(335*mm,433*mm,(-58.3+38.3)*mm/2.));       
402           fr1->AddNode(ppf,2,new TGeoTranslation(-335*mm,0.,(-58.3+38.3)*mm/2.));      fr1->AddNode(ppf,3,new TGeoTranslation(335*mm,0.,(-58.3+38.3)*mm/2.));
403           fr1->AddNode(ppf,0,new TGeoTranslation(-335*mm,-433*mm,(-58.3+38.3)*mm/2.)); fr1->AddNode(ppf,1,new TGeoTranslation(335*mm,-433*mm,(-58.3+38.3)*mm/2.));
404         
405           Double_t offsetx = 16.*mm, offsety = 34.*mm/2., interdistx = 48*mm+offsetx+0.6666*mm,interdisty = 19.*mm+2.*offsety;
406  
407        //gassiplex implementation
408        //it is in 3 different volumes: Hrec (in Hcufoil)+Hext   
409      
410           TGeoVolume *gassipl2 = gGeoManager->MakeBox("Hgassipl2",csi,32.*mm/2,3.*mm/2.,1.*mm/2.); //in Hrect
411           TGeoVolume *gassipl3 = gGeoManager->MakeBox("Hgassipl3",csi,60.*mm/2,3.*mm/2.,19.*mm/2.); //in Hfr1upcard
412           TGeoVolume *gassipl4 = gGeoManager->MakeBox("Hgassipl4",csi,60.*mm/2,3.*mm/2.,91.*mm/2.); //in Hext (the big rectangle of the card is 110 mm long, 62 mm wide and 1.5 mm high)
413           TGeoVolume *busext   = gGeoManager->MakeTubs("Hbusext",csi,29*mm,30*mm,40*mm/2.,0.,180); //in Hext
414           TGeoVolume *ext = new TGeoVolumeAssembly("Hext");
415           
416           rect->AddNode(gassipl2,1,new TGeoTranslation(0.,0.,0));
417           
418           for(Int_t hor=0; hor< 10; hor++){
419             for(Int_t vert=0; vert < 8; vert++){
420              cufoil->AddNode(rect,hor+vert*10,new TGeoTranslation(offsetx+ 48.*mm/2 + hor*interdistx-662.*mm/2,offsety + 19.*mm/2 + vert*interdisty-425.*mm/2.,0.));
421              fr1upcard->AddNode(gassipl3,hor+vert*10,new TGeoTranslation(offsetx+ 48.*mm/2 + hor*interdistx-662.*mm/2,offsety + 19.*mm/2 + vert*interdisty-425.*mm/2.,0.));
422              ext->AddNode(gassipl4,hor+vert*10,new TGeoTranslation(offsetx+ 48.*mm/2 + hor*interdistx-662.*mm/2,offsety + 19.*mm/2 +
423              vert*interdisty-425.*mm/2.,0));
424              ext->AddNode(busext,hor+vert*10,new TGeoTranslation(offsetx+ 48.*mm/2 + hor*interdistx-662.*mm/2,offsety + 19.*mm/2 +
425              vert*interdisty-425.*mm/2 + 3.*mm/2.,0));
426              }
427            }
428            
429            fr1up->AddNode(cufoil,4,new TGeoTranslation(-335*mm,433*mm,-20.0*mm/2+1.*mm/2));  fr1up->AddNode(cufoil,5,new TGeoTranslation(335*mm,433*mm,-20.0*mm/2+1.*mm/2));
430            fr1up->AddNode(cufoil,2,new TGeoTranslation(-335*mm,0,-20.0*mm/2+1.*mm/2));        fr1up->AddNode(cufoil,3,new TGeoTranslation(335*mm,0,-20.0*mm/2+1.*mm/2));
431            fr1up->AddNode(cufoil,0,new TGeoTranslation(-335*mm,-433*mm,-20.0*mm/2+1.*mm/2));  fr1up->AddNode(cufoil,1,new TGeoTranslation(335*mm,-433*mm,-20.0*mm/2+1.*mm/2));
432            
433            fr1up->AddNode(fr1upcard,4,new TGeoTranslation(-335*mm,433*mm,1.*mm/2.));   fr1up->AddNode(fr1upcard,5,new TGeoTranslation(335*mm,433*mm,1.*mm/2.));
434            fr1up->AddNode(fr1upcard,2,new TGeoTranslation(-335*mm,0,1.*mm/2.));        fr1up->AddNode(fr1upcard,3,new TGeoTranslation(335*mm,0,1.*mm/2.));
435            fr1up->AddNode(fr1upcard,0,new TGeoTranslation(-335*mm,-433*mm,1.*mm/2));  fr1up->AddNode(fr1upcard,1,new TGeoTranslation(335*mm,-433*mm,1.*mm/2.));
436
437
438   hmp->AddNode(ext,4,new TGeoTranslation(-335*mm,+433*mm,  (80.+1.7)*mm+58.3*mm+91*mm/2.)); hmp->AddNode(ext,5,new TGeoTranslation(+335*mm,+433*mm,  (80.+1.7)*mm+58.3*mm+91*mm/2.));
439   hmp->AddNode(ext,2,new TGeoTranslation(-335*mm,   0*mm,  (80.+1.7)*mm+58.3*mm+91*mm/2.)); hmp->AddNode(ext,3,new TGeoTranslation(+335*mm,   0*mm,  (80.+1.7)*mm+58.3*mm+91*mm/2.));
440   hmp->AddNode(ext,0,new TGeoTranslation(-335*mm,-433*mm,  (80.+1.7)*mm+58.3*mm+91*mm/2.)); hmp->AddNode(ext,1,new TGeoTranslation(+335*mm,-433*mm,  (80.+1.7)*mm+58.3*mm+91*mm/2.));
441
442         
443  hmp->AddNode(proxgap1,0,new TGeoTranslation(0.,0.,(9.-7.5)*mm/2.));//due to the TGeoVolumeAssembly definition the ch4 volume must be inserted around the collecting wires              
444  hmp->AddNode(proxgap2,0,new TGeoTranslation(0.,0.,(9+7.5 +34)*mm + (81.7-6.2-34.-9.-7.5)*mm/2.));// tot height(81.7) - Hsec - proxygap2 - top edge fr4 at (9+7.5) mm           
445                 
446 // ^ Y  single cell                                                5.5mm CH4 = 1*mm CsI + 4.45*mm CsI x cath +0.05*mm safety margin         
447 // |      ______________________________           
448 // |     |                              |          ^                            ||     
449 //       |                              | 1.05mm                                ||     
450 // 2.2*mm| xxxxxxxxxxxxxxxxxxxxxxxxxxxx |--              50um  x                || cat shift  x=0mm , y= 3.15mm , z=-2.70mm       
451 //       |                              |                                       ||     
452 //       |                              |                                       ||     
453 // __    |  ..........................  | 2.1mm                    20un .       ||  ano shift x=0mm , y= 2.00mm , z=-0.29mm   
454 //       |                              |                                       ||     
455 //       |                              |                                       ||     
456 //       | xxxxxxxxxxxxxxxxxxxxxxxxxxxx |--                    x                ||  cat shift x=0mm , y= 1.05mm , z=-2.70mm   
457 //       |                              |                                       ||     
458 //       |                              |         8.4mm                         ||   
459 // 4*mm  |                              | 2.1mm                                 ||  pad shift x=0mm , y= 0.00mm , z=2.25*mm   
460 //       |                              |                                       ||  
461 //       |                              |                                       ||  
462 //       | xxxxxxxxxxxxxxxxxxxxxxxxxxxx |--                    x                ||  cat shift x=0mm , y=-1.05mm , z=-2.70mm   
463 //       |                              |                                       ||  
464 //       |                              |                                       ||    
465 // __    |  ..........................  | 2.1mm                         . 2.04mm||  ano shift x=0mm , y=-2.00mm , z=-0.29mm   
466 //       |                              |                                       ||  
467 //       |                              |                                       ||  
468 //       | xxxxxxxxxxxxxxxxxxxxxxxxxxxx |--                    x    4.45mm      ||  cat shift x=0mm , y=-3.15mm , z=-2.70mm   
469 // 2.2*mm|                              |                                       ||  
470 //       |                              | 1.05mm                                ||         
471 //       |______________________________|          v                            ||    
472 //       <             8 mm             >                          
473 //                                   ----->X                                 ----->Z
474
475
476
477   ppf->AddNode(lar,0,new TGeoTranslation(-224.5*mm,-151.875*mm,  0.*mm));
478   ppf->AddNode(lar,1,new TGeoTranslation(-224.5*mm,- 50.625*mm,  0.*mm));
479   ppf->AddNode(lar,2,new TGeoTranslation(-224.5*mm,+ 50.625*mm,  0.*mm));
480   ppf->AddNode(lar,3,new TGeoTranslation(-224.5*mm,+151.875*mm,  0.*mm));
481   ppf->AddNode(lar,4,new TGeoTranslation(+224.5*mm,-151.875*mm,  0.*mm));
482   ppf->AddNode(lar,5,new TGeoTranslation(+224.5*mm,- 50.625*mm,  0.*mm));
483   ppf->AddNode(lar,6,new TGeoTranslation(+224.5*mm,+ 50.625*mm,  0.*mm));
484   ppf->AddNode(lar,7,new TGeoTranslation(+224.5*mm,+151.875*mm,  0.*mm));
485   ppf->AddNode(smo,0,new TGeoTranslation(- 65.0*mm,-151.875*mm,  0.*mm));
486   ppf->AddNode(smo,1,new TGeoTranslation(- 65.0*mm,- 50.625*mm,  0.*mm));
487   ppf->AddNode(smo,2,new TGeoTranslation(- 65.0*mm,+ 50.625*mm,  0.*mm));
488   ppf->AddNode(smo,3,new TGeoTranslation(- 65.0*mm,+151.875*mm,  0.*mm));
489   ppf->AddNode(smo,4,new TGeoTranslation(+ 65.0*mm,-151.875*mm,  0.*mm));
490   ppf->AddNode(smo,5,new TGeoTranslation(+ 65.0*mm,- 50.625*mm,  0.*mm));
491   ppf->AddNode(smo,6,new TGeoTranslation(+ 65.0*mm,+ 50.625*mm,  0.*mm));
492   ppf->AddNode(smo,7,new TGeoTranslation(+ 65.0*mm,+151.875*mm,  0.*mm)); 
493   
494
495 //hmp->AddNode(fr3,1,new TGeoTranslation(0.,0.,(81.7-29.)*mm-34.*mm/2));
496  hmp->AddNode(fr3,1,new TGeoTranslation(0.,0.,(9.+7.5)*mm+34.*mm/2));
497          fr3->AddNode( fr3up,1,    new TGeoTranslation(0.,  0.,  7*mm));
498          fr3->AddNode(fr3down,1,new TGeoTranslation(0.,  0., -10*mm));  
499
500 return hmp;
501   
502 }//CreateChamber()
503 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
504 void AliHMPIDv3::Init()
505 {
506 // This method defines ID for sensitive volumes, i.e. such geometry volumes for which there are if(TVirtualMC::GetMC()->CurrentVolID()==XXX) 
507 // statements in StepManager()
508 // Arguments: none
509 //   Returns: none      
510   AliDebug(1,"Start v2 HMPID.");    
511   //InitProperties();
512   AliDebug(1,"Stop v2 HMPID.");    
513 }
514 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
515 void AliHMPIDv3::DefineOpticalProperties() 
516 {
517   AliDebug(1,"");    
518
519 // Optical properties definition.
520   const Int_t kNbins=30;                 //number of photon energy points
521   Float_t emin=5.5,emax=8.5;             //Photon energy range,[eV]
522   Float_t deltaE = (emax - emin)/kNbins;
523   Float_t aEckov [kNbins]; 
524   Double_t dEckov [kNbins]; 
525   Float_t aAbsRad[kNbins], aAbsWin[kNbins], aAbsGap[kNbins], aAbsMet[kNbins];
526   Float_t aIdxRad[kNbins], aIdxWin[kNbins], aIdxGap[kNbins], aIdxMet[kNbins], aIdxPc[kNbins]; 
527   Float_t                                                    aQeAll [kNbins], aQePc [kNbins];
528   Double_t dReflMet[kNbins], dQePc[kNbins];
529
530   TF2 *pRaIF=new TF2("HidxRad","sqrt(1+0.554*(1239.84/x)^2/((1239.84/x)^2-5769)-0.0005*(y-20))"                                       ,emin,emax,0,50); //DiMauro mail temp 0-50 degrees C
531   TF1 *pWiIF=new TF1("HidxWin","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
532   TF1 *pGaIF=new TF1("HidxGap","1+0.12489e-6/(2.62e-4 - x*x/1239.84/1239.84)"                                                         ,emin,emax);      //?????? from where  
533
534   TF1 *pRaAF=new TF1("HabsRad","(x<7.8)*(gaus+gaus(3))+(x>=7.8)*0.0001"                                                               ,emin,emax);  //fit from DiMauro data 28.10.03 
535   pRaAF->SetParameters(3.20491e16,-0.00917890,0.742402,3035.37,4.81171,0.626309);
536   TF1 *pWiAF=new TF1("HabsWin","(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 
537   TF1 *pGaAF=new TF1("HabsGap","(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  
538   
539   TF1 *pQeF =new TF1("Hqe"    ,"0+(x>6.07267)*0.344811*(1-exp(-1.29730*(x-6.07267)))"                                                 ,emin,emax);  //fit from DiMauro data 28.10.03  
540                    
541   TString title=GetTitle();
542   Bool_t isFlatIdx=title.Contains("FlatIdx"); 
543   
544   for(Int_t i=0;i<kNbins;i++){
545     Float_t eV=emin+deltaE*i;  //Ckov energy in eV
546     aEckov [i] =1e-9*eV;    //Ckov energy in GeV
547     dEckov [i] = aEckov[i];
548     aAbsRad[i]=pRaAF->Eval(eV); (isFlatIdx)? aIdxRad[i]=1.292: aIdxRad[i]=pRaIF->Eval(eV,20);     
549     aAbsWin[i]=pWiAF->Eval(eV);              aIdxWin[i]=pWiIF->Eval(eV);
550     aAbsGap[i]=pGaAF->Eval(eV);              aIdxGap[i]=pGaIF->Eval(eV);   
551     aQeAll[i] =1;                     //QE for all other materials except for PC must be 1.  
552     aAbsMet[i] =0.0001;                aIdxMet[i]=0;                                             //metal ref idx must be 0 in order to reflect photon
553                                        aIdxPc [i]=1;           aQePc [i]=pQeF->Eval(eV);         //PC ref idx must be 1 in order to apply photon to QE conversion 
554     dQePc [i]= pQeF->Eval(eV);
555     dReflMet[i] = 0.;     // no reflection on the surface of the pc (?)                                       
556   }
557   TVirtualMC::GetMC()->SetCerenkov((*fIdtmed)[kC6F14]    , kNbins, aEckov, aAbsRad  , aQeAll , aIdxRad );    
558   TVirtualMC::GetMC()->SetCerenkov((*fIdtmed)[kSiO2]     , kNbins, aEckov, aAbsWin  , aQeAll , aIdxWin );    
559   TVirtualMC::GetMC()->SetCerenkov((*fIdtmed)[kCH4]      , kNbins, aEckov, aAbsGap  , aQeAll , aIdxGap );    
560   TVirtualMC::GetMC()->SetCerenkov((*fIdtmed)[kCu]       , kNbins, aEckov, aAbsMet  , aQeAll , aIdxMet );    
561   TVirtualMC::GetMC()->SetCerenkov((*fIdtmed)[kW]        , kNbins, aEckov, aAbsMet  , aQeAll , aIdxMet ); //n=0 means reflect photons       
562   TVirtualMC::GetMC()->SetCerenkov((*fIdtmed)[kCsI]      , kNbins, aEckov, aAbsMet  , aQePc  , aIdxPc  ); //n=1 means convert photons    
563   TVirtualMC::GetMC()->SetCerenkov((*fIdtmed)[kAl]       , kNbins, aEckov, aAbsMet  , aQeAll , aIdxMet );    
564
565   // Define a skin surface for the photocatode to enable 'detection' in G4
566   for(Int_t i=0; i<7; i++){
567   TVirtualMC::GetMC()->DefineOpSurface(Form("surfPc%i",i), kGlisur /*kUnified*/,kDielectric_metal,kPolished, 0.);
568   TVirtualMC::GetMC()->SetMaterialProperty(Form("surfPc%i",i), "EFFICIENCY", kNbins, dEckov, dQePc);
569   TVirtualMC::GetMC()->SetMaterialProperty(Form("surfPc%i",i), "REFLECTIVITY", kNbins, dEckov, dReflMet);
570   TVirtualMC::GetMC()->SetSkinSurface(Form("skinPc%i",i), Form("Hpad%i",i),Form("surfPc%i",i)); }
571
572   delete pRaAF;delete pWiAF;delete pGaAF; delete pRaIF; delete pWiIF; delete pGaIF; delete pQeF;
573 }
574 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
575 Bool_t AliHMPIDv3::IsLostByFresnel()
576 {
577 // Calculate probability for the photon to be lost by Fresnel reflection.
578   TLorentzVector p4;
579   Double_t mom[3],localMom[3];
580   TVirtualMC::GetMC()->TrackMomentum(p4);   mom[0]=p4(1);   mom[1]=p4(2);   mom[2]=p4(3);
581   localMom[0]=0; localMom[1]=0; localMom[2]=0;
582   TVirtualMC::GetMC()->Gmtod(mom,localMom,2);
583   Double_t localTc    = localMom[0]*localMom[0]+localMom[2]*localMom[2];
584   Double_t localTheta = TMath::ATan2(TMath::Sqrt(localTc),localMom[1]);
585   Double_t cotheta = TMath::Abs(TMath::Cos(localTheta));
586   if(TVirtualMC::GetMC()->GetRandom()->Rndm() < Fresnel(p4.E()*1e9,cotheta,1)){
587     AliDebug(1,"Photon lost");
588     return kTRUE;
589   }else
590     return kFALSE;
591 }//IsLostByFresnel()
592 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
593 void AliHMPIDv3::GenFee(Float_t qtot)
594 {
595 // Generate FeedBack photons for the current particle. To be invoked from StepManager().
596 // eloss=0 means photon so only pulse height distribution is to be analysed.
597   TLorentzVector x4;
598   TVirtualMC::GetMC()->TrackPosition(x4); 
599   Int_t iNphotons=TVirtualMC::GetMC()->GetRandom()->Poisson(0.02*qtot);  //# of feedback photons is proportional to the charge of hit
600   AliDebug(1,Form("N photons=%i",iNphotons));
601   Int_t j;
602   Float_t cthf, phif, enfp = 0, sthf, e1[3], e2[3], e3[3], vmod, uswop,dir[3], phi,pol[3], mom[4];
603 //Generate photons
604   for(Int_t i=0;i<iNphotons;i++){//feedbacks loop
605     Double_t ranf[2];
606     TVirtualMC::GetMC()->GetRandom()->RndmArray(2,ranf);    //Sample direction
607     cthf=ranf[0]*2-1.0;
608     if(cthf<0) continue;
609     sthf = TMath::Sqrt((1. - cthf) * (1. + cthf));
610     phif = ranf[1] * 2 * TMath::Pi();
611     
612     if(Double_t randomNumber=TVirtualMC::GetMC()->GetRandom()->Rndm()<=0.57)
613       enfp = 7.5e-9;
614     else if(randomNumber<=0.7)
615       enfp = 6.4e-9;
616     else
617       enfp = 7.9e-9;
618     
619
620     dir[0] = sthf * TMath::Sin(phif);    dir[1] = cthf;    dir[2] = sthf * TMath::Cos(phif);
621     TVirtualMC::GetMC()->Gdtom(dir, mom, 2);
622     mom[0]*=enfp;    mom[1]*=enfp;    mom[2]*=enfp;
623     mom[3] = TMath::Sqrt(mom[0]*mom[0]+mom[1]*mom[1]+mom[2]*mom[2]);
624     
625     // Polarisation
626     e1[0]=      0;    e1[1]=-dir[2];    e1[2]= dir[1];
627     e2[0]=-dir[1];    e2[1]= dir[0];    e2[2]=      0;
628     e3[0]= dir[1];    e3[1]=      0;    e3[2]=-dir[0];
629     
630     vmod=0;
631     for(j=0;j<3;j++) vmod+=e1[j]*e1[j];
632     if (!vmod) for(j=0;j<3;j++) {
633       uswop=e1[j];
634       e1[j]=e3[j];
635       e3[j]=uswop;
636     }
637     vmod=0;
638     for(j=0;j<3;j++) vmod+=e2[j]*e2[j];
639     if (!vmod) for(j=0;j<3;j++) {
640       uswop=e2[j];
641       e2[j]=e3[j];
642       e3[j]=uswop;
643     }
644     
645     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;    
646     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;
647     
648     phi = TVirtualMC::GetMC()->GetRandom()->Rndm()* 2 * TMath::Pi();
649     for(j=0;j<3;j++) pol[j]=e1[j]*TMath::Sin(phi)+e2[j]*TMath::Cos(phi);
650     TVirtualMC::GetMC()->Gdtom(pol, pol, 2);
651     Int_t outputNtracksStored;    
652     gAlice->GetMCApp()->PushTrack(1,                             //transport
653                      gAlice->GetMCApp()->GetCurrentTrackNumber(),//parent track 
654                      50000051,                                   //PID
655                      mom[0],mom[1],mom[2],mom[3],                //track momentum  
656                      x4.X(),x4.Y(),x4.Z(),x4.T(),                //track origin 
657                      pol[0],pol[1],pol[2],                       //polarization
658                      kPFeedBackPhoton,                           //process ID   
659                      outputNtracksStored,                        //on return how many new photons stored on stack
660                      1.0);                                       //weight
661   }//feedbacks loop
662   AliDebug(1,"Stop.");
663 }//GenerateFeedbacks()
664 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
665 void AliHMPIDv3::Hits2SDigits()
666 {
667 // Interface method ivoked from AliSimulation to create a list of sdigits corresponding to list of hits. Every hit generates one or more sdigits.
668 // Arguments: none
669 //   Returns: none   
670   AliDebug(1,"Start.");
671   for(Int_t iEvt=0;iEvt < GetLoader()->GetRunLoader()->GetNumberOfEvents();iEvt++){                //events loop
672     GetLoader()->GetRunLoader()->GetEvent(iEvt);                          //get next event
673   
674     if(!GetLoader()->TreeH()) {GetLoader()->LoadHits();                    }
675     if(!GetLoader()->TreeS()) {GetLoader()->MakeTree("S"); MakeBranch("S");}//to
676           
677     for(Int_t iEnt=0;iEnt<GetLoader()->TreeH()->GetEntries();iEnt++){//prims loop
678       GetLoader()->TreeH()->GetEntry(iEnt);
679       Hit2Sdi(Hits(),SdiLst());
680     }//prims loop
681     GetLoader()->TreeS()->Fill();
682     GetLoader()->WriteSDigits("OVERWRITE");
683     SdiReset();
684   }//events loop  
685   GetLoader()->UnloadHits();
686   GetLoader()->UnloadSDigits();  
687   AliDebug(1,"Stop.");
688 }//Hits2SDigits()
689 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
690 void AliHMPIDv3::Hit2Sdi(TClonesArray *pHitLst,TClonesArray *pSdiLst)
691 {
692 // Converts list of hits to list of sdigits. 
693 // Arguments: pHitLst  - list of hits provided not empty
694 //            pSDigLst - list of sdigits where to store the results
695 //   Returns: none         
696   for(Int_t iHit=0;iHit<pHitLst->GetEntries();iHit++){         //hits loop
697     AliHMPIDHit *pHit=(AliHMPIDHit*)pHitLst->At(iHit);         //get pointer to current hit   
698     pHit->Hit2Sdi(pSdiLst);                                    //convert this hit to list of sdigits     
699   }//hits loop loop
700 }//Hits2Sdi()
701 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
702 void AliHMPIDv3::Digits2Raw()
703 {
704 // Interface method invoked by AliSimulation to create raw data streams from digits. Events loop is done in AliSimulation
705 // Arguments: none
706 //   Returns: none    
707   AliDebug(1,"Start.");
708   GetLoader()->LoadDigits();
709   TTree * treeD = GetLoader()->TreeD();
710   if(!treeD) {
711     AliError("No digits tree!");
712     return;
713   }
714   treeD->GetEntry(0);
715   
716   
717   AliHMPIDRawStream *pRS=0x0;
718   pRS->WriteRaw(DigLst());
719    
720   GetLoader()->UnloadDigits();
721   AliDebug(1,"Stop.");      
722 }//Digits2Raw()
723 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
724 Float_t AliHMPIDv3::Fresnel(Float_t ene,Float_t pdoti, Bool_t pola)
725 {
726 // Correction for Fresnel   ???????????
727 // Arguments:   ene - photon energy [GeV],
728 //              PDOTI=COS(INC.ANG.), PDOTR=COS(POL.PLANE ROT.ANG.)
729 //   Returns:  
730     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,
731                       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,
732                       7.8,7.9,8.0,8.1,8.2,8.3,8.4,8.5};
733     Float_t csin[36] = {2.14,2.21,2.33,2.48,2.76,2.97,2.99,2.59,2.81,3.05,
734                         2.86,2.53,2.55,2.66,2.79,2.96,3.18,3.05,2.84,2.81,2.38,2.11,
735                         2.01,2.13,2.39,2.73,3.08,3.15,2.95,2.73,2.56,2.41,2.12,1.95,
736                         1.72,1.53};
737     Float_t csik[36] = {0.,0.,0.,0.,0.,0.196,0.408,0.208,0.118,0.49,0.784,0.543,
738                         0.424,0.404,0.371,0.514,0.922,1.102,1.139,1.376,1.461,1.253,0.878,
739                         0.69,0.612,0.649,0.824,1.347,1.571,1.678,1.763,1.857,1.824,1.824,
740                         1.714,1.498};
741     Float_t xe=ene;
742     Int_t  j=Int_t(xe*10)-49;
743     Float_t cn=csin[j]+((csin[j+1]-csin[j])/0.1)*(xe-en[j]);
744     Float_t ck=csik[j]+((csik[j+1]-csik[j])/0.1)*(xe-en[j]);
745
746     //FORMULAE FROM HANDBOOK OF OPTICS, 33.23 OR
747     //W.R. HUNTER, J.O.S.A. 54 (1964),15 , J.O.S.A. 55(1965),1197
748
749     Float_t sinin=TMath::Sqrt((1.-pdoti)*(1.+pdoti));
750     Float_t tanin=sinin/pdoti;
751
752     Float_t c1=cn*cn-ck*ck-sinin*sinin;
753     Float_t c2=4*cn*cn*ck*ck;
754     Float_t aO=TMath::Sqrt(0.5*(TMath::Sqrt(c1*c1+c2)+c1));
755     Float_t b2=0.5*(TMath::Sqrt(c1*c1+c2)-c1);
756     
757     Float_t rs=((aO-pdoti)*(aO-pdoti)+b2)/((aO+pdoti)*(aO+pdoti)+b2);
758     Float_t rp=rs*((aO-sinin*tanin)*(aO-sinin*tanin)+b2)/((aO+sinin*tanin)*(aO+sinin*tanin)+b2);
759     
760
761     //CORRECTION FACTOR FOR SURFACE ROUGHNESS
762     //B.J. STAGG  APPLIED OPTICS, 30(1991),4113
763
764     Float_t sigraf=18.;
765     Float_t lamb=1240/ene;
766     Float_t fresn;
767  
768     Float_t  rO=TMath::Exp(-(4*TMath::Pi()*pdoti*sigraf/lamb)*(4*TMath::Pi()*pdoti*sigraf/lamb));
769
770     if(pola)
771     {
772         Float_t pdotr=0.8;                                 //DEGREE OF POLARIZATION : 1->P , -1->S
773         fresn=0.5*(rp*(1+pdotr)+rs*(1-pdotr));
774     }
775     else
776         fresn=0.5*(rp+rs);
777       
778     fresn = fresn*rO;
779     return fresn;
780 }//Fresnel()
781 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
782 void AliHMPIDv3::Print(Option_t *option)const
783 {
784 // Debug printout
785   TObject::Print(option);
786 }//void AliHMPID::Print(Option_t *option)const
787 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
788
789 Bool_t AliHMPIDv3::Raw2SDigits(AliRawReader *pRR)
790 {
791 // Arguments: pRR- raw reader 
792 //   Returns: kTRUE on success (currently ignored in AliSimulation::ConvertRaw2SDigits())      
793 //AliHMPIDDigit sdi; //tmp sdigit, raw digit will be converted to it
794   
795   if(!GetLoader()->TreeS()) {MakeTree("S");  MakeBranch("S");}
796     
797   TClonesArray *pSdiLst=SdiLst(); Int_t iSdiCnt=0; //tmp list of sdigits for all chambers
798   AliHMPIDRawStream stream(pRR);
799   while(stream.Next())
800   {
801     for(Int_t iPad=0;iPad<stream.GetNPads();iPad++) {
802       AliHMPIDDigit sdi(stream.GetPadArray()[iPad],stream.GetChargeArray()[iPad]);
803       new((*pSdiLst)[iSdiCnt++]) AliHMPIDDigit(sdi); //add this digit to the tmp list
804     }
805   }
806   
807   GetLoader()->TreeS()->Fill(); GetLoader()->WriteSDigits("OVERWRITE");//write out sdigits
808   SdiReset();
809   return kTRUE;
810
811 }//Raw2SDigits
812
813 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
814 void AliHMPIDv3::StepCount()
815 {
816 // Count number of ckovs created  
817 }
818 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
819 void AliHMPIDv3::StepHistory()
820 {
821 // This methode is invoked from StepManager() in order to print out 
822   static Int_t iStepN;
823   const char *sParticle;
824   switch(TVirtualMC::GetMC()->TrackPid()){
825     case kProton:      sParticle="PROTON"    ;break;
826     case kNeutron:     sParticle="neutron"   ;break;
827     case kGamma:       sParticle="gamma"     ;break;
828     case 50000050:     sParticle="CKOV"      ;break;
829     case kPi0:         sParticle="Pi0"       ;break;  
830     case kPiPlus:      sParticle="Pi+"       ;break;  
831     case kPiMinus:     sParticle="Pi-"       ;break;  
832     case kElectron:    sParticle="electron"  ;break;  
833     default:           sParticle="not known" ;break;
834   }
835
836   TString flag="fanny combination";
837   if(TVirtualMC::GetMC()->IsTrackAlive()) {
838     if(TVirtualMC::GetMC()->IsTrackEntering())      flag="enters to";
839     else if(TVirtualMC::GetMC()->IsTrackExiting())  flag="exits from";
840     else if(TVirtualMC::GetMC()->IsTrackInside())   flag="inside";
841   } else {
842     if(TVirtualMC::GetMC()->IsTrackStop())          flag="stopped in";
843   }
844   
845   Int_t vid=0,copy=0;
846   TString path=TVirtualMC::GetMC()->CurrentVolName(); path.Prepend("-");path.Prepend(TVirtualMC::GetMC()->CurrentVolOffName(1));//current volume and his mother are always there
847   vid=TVirtualMC::GetMC()->CurrentVolOffID(2,copy);  if(vid) {path.Prepend("-");path.Prepend(TVirtualMC::GetMC()->VolName(vid));}
848   vid=TVirtualMC::GetMC()->CurrentVolOffID(3,copy);  if(vid) {path.Prepend("-");path.Prepend(TVirtualMC::GetMC()->VolName(vid));}
849  
850   
851   Printf("Step %i: %s (%i) %s %s m=%.6f GeV q=%.1f dEdX=%.4f Etot=%.4f",iStepN,sParticle,TVirtualMC::GetMC()->TrackPid(),flag.Data(),path.Data(),TVirtualMC::GetMC()->TrackMass(),TVirtualMC::GetMC()->TrackCharge(),TVirtualMC::GetMC()->Edep()*1e9,TVirtualMC::GetMC()->Etot());
852   
853   Double_t gMcTrackPos[3]; TVirtualMC::GetMC()->TrackPosition(gMcTrackPos[0],gMcTrackPos[1],gMcTrackPos[2]);
854   Double_t  gMcTrackPosLoc[3]; TVirtualMC::GetMC()->Gmtod(gMcTrackPos,gMcTrackPosLoc,1);
855   Printf("TVirtualMC::GetMC() Track Position (MARS) x: %5.3lf, y: %5.3lf, z: %5.3lf (r: %5.3lf) ---> (LOC) x: %5.3f, y: %5.3f, z: %5.3f",gMcTrackPos[0],gMcTrackPos[1],gMcTrackPos[2],TMath::Sqrt(gMcTrackPos[0]*gMcTrackPos[0]+gMcTrackPos[1]*gMcTrackPos[1]+gMcTrackPos[2]*gMcTrackPos[2]),gMcTrackPosLoc[0],gMcTrackPosLoc[1],gMcTrackPosLoc[2]);
856   
857
858   
859   Printf("Step %i: tid=%i flags alive=%i disap=%i enter=%i exit=%i inside=%i out=%i stop=%i new=%i",
860                             iStepN, gAlice->GetMCApp()->GetCurrentTrackNumber(),
861                             TVirtualMC::GetMC()->IsTrackAlive(), TVirtualMC::GetMC()->IsTrackDisappeared(),TVirtualMC::GetMC()->IsTrackEntering(), TVirtualMC::GetMC()->IsTrackExiting(),
862                             TVirtualMC::GetMC()->IsTrackInside(),TVirtualMC::GetMC()->IsTrackOut(),        TVirtualMC::GetMC()->IsTrackStop(),     TVirtualMC::GetMC()->IsNewTrack());
863   
864   Float_t a,z,den,rad,abs; a=z=den=rad=abs=-1;
865   Int_t mid=TVirtualMC::GetMC()->CurrentMaterial(a,z,den,rad,abs);
866   Printf("Step %i: mid=%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);
867   
868   TArrayI proc;  TVirtualMC::GetMC()->StepProcesses(proc); 
869   Printf("Processes in this step:");
870   for ( int i = 0 ; i < proc.GetSize(); i++)
871   {
872     Printf("%s",TMCProcessName[proc.At(i)]);
873   }
874   Printf("End process list");
875   
876   iStepN++;
877 }//StepHistory()
878 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
879 void AliHMPIDv3::StepManager()
880 {
881 // Full Step Manager.
882 // Arguments: none
883 //   Returns: none           
884 //  StepHistory(); return; //uncomment to print tracks history
885  //  StepCount(); return;     //uncomment to count photons
886   
887    TString volname = TVirtualMC::GetMC()->CurrentVolName();
888
889 //Treat photons    
890     if((TVirtualMC::GetMC()->TrackPid()==50000050||TVirtualMC::GetMC()->TrackPid()==50000051)&&volname.Contains("Hpad")){ //photon (Ckov or feedback) hits on module PC (Hpad)
891     if(TVirtualMC::GetMC()->Edep()>0){                                                                           //photon survided QE test i.e. produces electron
892       if(IsLostByFresnel()){ TVirtualMC::GetMC()->StopTrack(); return;}                                          //photon lost due to fersnel reflection on PC       
893       Int_t   tid=     TVirtualMC::GetMC()->GetStack()->GetCurrentTrackNumber();                                 //take TID
894       Int_t   pid=     TVirtualMC::GetMC()->TrackPid();                                                          //take PID
895       Float_t etot=    TVirtualMC::GetMC()->Etot();                                                              //total hpoton energy, [GeV] 
896       Double_t x[3];   TVirtualMC::GetMC()->TrackPosition(x[0],x[1],x[2]);                                       //take MARS position at entrance to PC
897       Float_t hitTime= (Float_t)TVirtualMC::GetMC()->TrackTime();                                                //hit formation time       
898       TString tmpname = volname; tmpname.Remove(0,4); Int_t idch = tmpname.Atoi();               //retrieve the chamber number
899       Float_t xl,yl;   AliHMPIDParam::Instance()->Mars2Lors(idch,x,xl,yl);                       //take LORS position 
900       new((*fHits)[fNhits++])AliHMPIDHit(idch,etot,pid,tid,xl,yl,hitTime,x);                             //HIT for photon, position at P, etot will be set to Q
901       if(fDoFeed) GenFee(etot);                                                                  //generate feedback photons etot is modified in hit ctor to Q of hit
902     }//photon hit PC and DE >0 
903   }//photon hit PC
904    
905   
906   //Treat charged particles  
907   static Float_t eloss;                                                                           //need to store mip parameters between different steps    
908   static Double_t in[3];                                                                          
909
910   if(TVirtualMC::GetMC()->IsTrackEntering() && TVirtualMC::GetMC()->TrackCharge() && volname.Contains("Hpad")) //Trackref stored when entering in the pad volume
911     AddTrackReference(TVirtualMC::GetMC()->GetStack()->GetCurrentTrackNumber(), AliTrackReference::kHMPID);       //for acceptance calculations
912    
913
914   if(TVirtualMC::GetMC()->TrackCharge() && volname.Contains("Hcel")){                                           //charged particle in amplification gap (Hcel)
915     if(TVirtualMC::GetMC()->IsTrackEntering()||TVirtualMC::GetMC()->IsNewTrack()) {                                               //entering or newly created
916       eloss=0;                                                                                    //reset Eloss collector                         
917       TVirtualMC::GetMC()->TrackPosition(in[0],in[1],in[2]);                                                      //take position at the entrance
918     }else if(TVirtualMC::GetMC()->IsTrackExiting()||TVirtualMC::GetMC()->IsTrackStop()||TVirtualMC::GetMC()->IsTrackDisappeared()){               //exiting or disappeared
919       eloss              +=TVirtualMC::GetMC()->Edep();                                                           //take into account last step Eloss
920       Int_t tid=          TVirtualMC::GetMC()->GetStack()->GetCurrentTrackNumber();                               //take TID
921       Int_t pid=          TVirtualMC::GetMC()->TrackPid();                                                        //take PID
922       Double_t out[3];    TVirtualMC::GetMC()->TrackPosition(out[0],out[1],out[2]);                               //take MARS position at exit
923       Float_t hitTime= (Float_t)TVirtualMC::GetMC()->TrackTime();                                                         //hit formation time       
924       out[0]=0.5*(out[0]+in[0]);                                                                  //
925       out[1]=0.5*(out[1]+in[1]);                                                                  //take hit position at the anod plane
926       out[2]=0.5*(out[2]+in[2]);
927       TString tmpname = volname;  tmpname.Remove(0,4);  Int_t idch = tmpname.Atoi();              //retrieve the chamber number
928       Float_t xl,yl;AliHMPIDParam::Instance()->Mars2Lors(idch,out,xl,yl);                         //take LORS position
929       if(eloss>0) {
930         new((*fHits)[fNhits++])AliHMPIDHit(idch,eloss,pid,tid,xl,yl,hitTime,out);                           //HIT for MIP, position near anod plane, eloss will be set to Q 
931         if(fDoFeed) GenFee(eloss);                                                                  //generate feedback photons 
932       }
933     }else                                                                                         //just going inside
934       eloss          += TVirtualMC::GetMC()->Edep();                                                              //collect this step eloss 
935   }//MIP in GAP
936  
937 }//StepManager()
938 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
939 void AliHMPIDv3::TestPoint(Int_t ch,Float_t x,Float_t y)
940 {
941 // Utility method to check the validity of geometry by poviding some crucial points
942 // Arguments: ch,x,y- crucial point definition (cm) in LORS
943 //   Returns: none    
944   Double_t mars[3];
945   AliHMPIDParam::Instance()->Lors2Mars(ch,x,y,mars);
946   Printf("(ch=%i,locX=%.2f,locY=%.2f) %s",ch,x,y,gGeoManager->FindNode(mars[0],mars[1],mars[2])->GetName());
947 }//TestPoint()
948 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
949 void AliHMPIDv3::TestGeom()
950 {
951 //
952 // Test method to check geometry
953 //
954   //TGeoManager::Import("misaligned_geometry.root");
955   TGeoManager::Import("geometry.root");
956   for(Int_t ch=AliHMPIDParam::kMinCh;ch<=AliHMPIDParam::kMaxCh;ch++)
957     TestPoint(ch,0,0);
958 }//TestPoint()
959 //+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
960
961 void  AliHMPIDv3::IdealPosition(Int_t iCh,TGeoHMatrix *pMatrix)       //ideal position of given chamber 
962 {
963 // Construct ideal position matrix for a given chamber
964 // Arguments: iCh- chamber ID; pMatrix- pointer to precreated unity matrix where to store the results
965 //   Returns: none
966   const Double_t kAngHor=19.5;        //  horizontal angle between chambers  19.5 grad
967   const Double_t kAngVer=20;          //  vertical angle between chambers    20   grad     
968   const Double_t kAngCom=30;          //  common HMPID rotation with respect to x axis  30   grad     
969   const Double_t kTrans[3]={490,0,0}; //  center of the chamber is on window-gap surface
970   pMatrix->RotateY(90);               //  rotate around y since initial position is in XY plane -> now in YZ plane
971   pMatrix->SetTranslation(kTrans);    //  now plane in YZ is shifted along x 
972   switch(iCh){
973     case 0:                pMatrix->RotateY(kAngHor);  pMatrix->RotateZ(-kAngVer);  break; //right and down 
974     case 1:                                            pMatrix->RotateZ(-kAngVer);  break; //down              
975     case 2:                pMatrix->RotateY(kAngHor);                               break; //right 
976     case 3:                                                                         break; //no rotation
977     case 4:                pMatrix->RotateY(-kAngHor);                              break; //left   
978     case 5:                                            pMatrix->RotateZ(kAngVer);   break; //up
979     case 6:                pMatrix->RotateY(-kAngHor); pMatrix->RotateZ(kAngVer);   break; //left and up 
980   }
981   pMatrix->RotateZ(kAngCom);     //apply common rotation  in XY plane    
982 }
983
984
985 void  AliHMPIDv3::IdealPositionCradle(Int_t iCh,TGeoHMatrix *pMatrix)     //ideal position of given one module of the cradle
986 {
987 // Construct ideal position matrix for a given module cradle
988 // Arguments: iCh- chamber ID; pMatrix- pointer to precreated unity matrix where to store the results
989 //   Returns: none
990   const Double_t kAngHor=19.5;        //  horizontal angle between chambers  19.5 grad
991   const Double_t kAngVer=20;          //  vertical angle between chambers    20   grad
992   const Double_t kAngCom=30;          //  common HMPID rotation with respect to x axis  30   grad
993   const Double_t kTrans[3]={423.+ 29,0,67}; //  z-center of the cradle module
994   pMatrix->RotateY(90);               //  rotate around y since initial position is in XY plane -> now in YZ plane
995   pMatrix->SetTranslation(kTrans);    //  now plane in YZ is shifted along x
996   switch(iCh){
997     case 0:                pMatrix->RotateY(kAngHor);  pMatrix->RotateZ(-kAngVer);  break; //right and down
998     case 1:                                            pMatrix->RotateZ(-kAngVer);  break; //down
999     case 2:                pMatrix->RotateY(kAngHor);                               break; //right
1000     case 3:                                                                         break; //no rotation
1001     case 4:                pMatrix->RotateY(-kAngHor);                              break; //left
1002     case 5:                                            pMatrix->RotateZ(kAngVer);   break; //up
1003     case 6:                pMatrix->RotateY(-kAngHor); pMatrix->RotateZ(kAngVer);   break; //left and up
1004   }
1005   pMatrix->RotateZ(kAngCom);     //apply common rotation  in XY plane
1006 }
1007
1008
1009
1010 TGeoVolume* AliHMPIDv3::CreateCradle()
1011 {
1012
1013 //Method that builds the Cradle geometry
1014 //according to the base topology created
1015 //in CradleBaseVolume(...)
1016
1017   Double_t mm = 0.1;
1018
1019   Double_t params[10]={0.5,10.,24.,-1,5.2,1.5,3.5,8.5,3.8,0.};
1020   TGeoMedium *med   =gGeoManager->GetMedium("HMPID_Al");  
1021   TGeoVolume *cradle=new TGeoVolumeAssembly("Hcradle");
1022   
1023   //Double_t baselong[7]={6037*mm-2*60*mm, 6037*mm-2*60*mm,60*mm,0.,100*mm,10*mm,10*mm};//2CRE2112P3
1024   Double_t baselong[7]={6037*mm-2*100*mm, 6037*mm-2*100*mm,60*mm,0.,100*mm,10*mm,10*mm};//2CRE2112P3
1025   TGeoVolume *lbase = CradleBaseVolume(med,baselong,"cradleLbase");
1026   lbase->SetLineColor(kGray);
1027
1028   Double_t baseshort[7]={1288.*mm+2*100*mm, 1288.*mm+2*100*mm,60*mm,0.,100*mm,10*mm,10*mm};//2CRE2112P3
1029   TGeoVolume *sbase = CradleBaseVolume(med,baseshort,"cradleSbase");
1030   sbase->SetLineColor(kGray);
1031
1032   //one side
1033
1034   Double_t height = 30.*mm; //30 = 2*(1488/2-729) (2CRE2112P3)
1035   Double_t tubeheight = 50.*mm; Double_t heightred = 5.*mm; Double_t zred = 5.*mm;
1036   Double_t oneshift = tubeheight/TMath::Tan(TMath::DegToRad()*20.)+(1458.-35)*mm/2 - (1607-35)*mm/2;
1037   Double_t linclined[7] = {1458.*mm-params[6]-0.5,1607.*mm-params[6]-0.5,tubeheight,oneshift, height ,heightred,zred}; //3.5 is for not correct measurements in 2CRE2112P3<=> 597!=inclined*sin(20) 
1038   TGeoVolume *inclin = CradleBaseVolume(med,linclined,"inclinedbar");
1039   inclin->SetLineColor(kGray);
1040   Double_t lhorizontal[7] = {1641.36*mm+params[7],1659.*mm+params[7],tubeheight,0, height ,heightred,zred};
1041   TGeoVolume *horiz = CradleBaseVolume(med,lhorizontal,"horizontalbar");
1042   horiz->SetLineColor(kGray);
1043   
1044   //inner bars, they are named as the numbering in 2CRE2112P3
1045   Double_t fourshift = tubeheight/TMath::Tan(TMath::DegToRad()*55.);  
1046   Double_t lfour[7] = {592*mm,592*mm,tubeheight,fourshift,height,heightred,zred};  
1047   TGeoVolume *four = CradleBaseVolume(med,lfour,"bar4");
1048   four->SetLineColor(kGray);
1049
1050   Double_t fiveshift = tubeheight/TMath::Tan(TMath::DegToRad()*75);
1051   Double_t lfive[7] = {500.*mm,500.*mm,tubeheight,fiveshift,height,heightred,zred};  
1052   TGeoVolume *five = CradleBaseVolume(med,lfive,"bar5");
1053   five->SetLineColor(kGray);
1054   
1055   Double_t sixshift = tubeheight/TMath::Tan(TMath::DegToRad()*55)+459*mm/2-480*mm/2;  
1056   Double_t lsix[7] = {456*mm,477*mm,tubeheight,sixshift,height,heightred,zred};  
1057   TGeoVolume *six = CradleBaseVolume(med,lsix,"bar6");
1058   six->SetLineColor(kGray);
1059   
1060   Double_t sevenshift = tubeheight/TMath::Tan(TMath::DegToRad()*50)+472*mm/2-429.*mm/2;
1061   Double_t lseven[7] = {429*mm,472*mm,tubeheight,sevenshift,height,heightred,zred};  
1062   TGeoVolume *seven = CradleBaseVolume(med,lseven,"bar7");
1063   seven->SetLineColor(kGray);
1064     
1065   Double_t eightshift = tubeheight/TMath::Tan(TMath::DegToRad()*30)+244.*mm/2-200.*mm/2 -3;
1066   Double_t leight[7] = {200.*mm,244.*mm,tubeheight,eightshift,height,heightred,zred};  
1067   TGeoVolume *eight = CradleBaseVolume(med,leight,"bar8");
1068   eight->SetLineColor(kGray);
1069
1070   Double_t nineshift = -tubeheight/TMath::Tan(TMath::DegToRad()*71)+83.*mm/2-66.*mm/2;
1071   Double_t lnine[7] = {59.5*mm,76.5*mm,tubeheight,nineshift,height,heightred,zred};  
1072   TGeoVolume *nine = CradleBaseVolume(med,lnine,"bar9");
1073   nine->SetLineColor(kGray);
1074
1075   Double_t tenshift = (-tubeheight/TMath::Tan(TMath::DegToRad()*60) -221.*mm/2+195.*mm/2);
1076   Double_t lten[7] = {195.*mm,221.*mm,tubeheight,tenshift,height,heightred,zred};  
1077   TGeoVolume *ten = CradleBaseVolume(med,lten,"bar10");
1078   ten->SetLineColor(kGray);
1079   
1080   Double_t elevenshift = (-tubeheight/TMath::Tan(TMath::DegToRad()*70) -338.*mm/2+315.*mm/2);
1081   Double_t leleven[7] = {308.*mm,331.*mm,tubeheight,elevenshift,height,heightred,zred};  
1082   TGeoVolume *eleven = CradleBaseVolume(med,leleven,"bar11");
1083   eleven->SetLineColor(kGray);
1084     
1085   Double_t twelveshift = (-tubeheight/TMath::Tan(TMath::DegToRad()*60) -538.*mm/2+508.*mm/2);
1086   Double_t ltwelve[7] = {507.*mm,537.*mm,tubeheight,twelveshift,height,heightred,zred};  
1087   TGeoVolume *twelve = CradleBaseVolume(med,ltwelve,"bar12");
1088   twelve->SetLineColor(kGray);
1089
1090   Double_t thirteenshift = tubeheight/TMath::Tan(TMath::DegToRad()*43); 
1091   Double_t lthirteen[7] = {708.*mm,708.*mm,tubeheight,thirteenshift,height,heightred,zred};  
1092   TGeoVolume *thirteen = CradleBaseVolume(med,lthirteen,"bar13");
1093   thirteen->SetLineColor(kGray); 
1094   
1095   
1096   //vertical rectangles
1097    TGeoVolume *vbox= new TGeoVolumeAssembly("Hvbox");
1098    vbox->SetLineColor(kViolet);
1099    Double_t width = 50.*mm;
1100    
1101    TGeoVolume *vboxlast= new TGeoVolumeAssembly("Hvboxlast");//vertical structure on the short base
1102    vboxlast->SetLineColor(kViolet);
1103   
1104   Double_t barheight = 100.*mm; 
1105   Double_t lAfourteen[7] = {1488.*mm,1488.*mm,barheight,0,width,heightred,zred};  
1106   TGeoVolume *afourteen = CradleBaseVolume(med,lAfourteen,"bar14top");
1107   afourteen->SetLineColor(kGray); 
1108  
1109   Double_t lBfourteen[7] = {387*mm,387.*mm,barheight,0,width,heightred,zred};  
1110   TGeoVolume *bfourteen = CradleBaseVolume(med,lBfourteen,"bar14vert");
1111   bfourteen->SetLineColor(kGray);
1112    
1113   Double_t lCfourteen[7] = {1288.*mm,1288.*mm,barheight,0,width,heightred,zred};  
1114   TGeoVolume *cfourteen = CradleBaseVolume(med,lCfourteen,"bar14bot");
1115   cfourteen->SetLineColor(kGray);
1116
1117   Double_t oblshift = 50.*mm/ TMath::Tan(TMath::DegToRad()*35); 
1118   Double_t lDfourteen[7] = {603.*mm,603.*mm,50.*mm,oblshift,width,heightred,zred}; 
1119   TGeoVolume *dfourteen = CradleBaseVolume(med,lDfourteen,"bar14incl");
1120   dfourteen->SetLineColor(kGray);
1121   
1122   
1123   Double_t lDfourteenlast[7] = {667.*mm,667.*mm,50.*mm,oblshift,width,heightred,zred};  
1124   TGeoVolume *dfourteenlast = CradleBaseVolume(med,lDfourteenlast,"bar14incllast");
1125   dfourteenlast->SetLineColor(kGray);
1126   
1127   vbox->AddNode(afourteen,1,new TGeoTranslation(0.,487.*mm/2 -100.*mm/2,0.));
1128   TGeoRotation *vinrot = new TGeoRotation("vertbar"); vinrot->RotateZ(90);
1129   vbox->AddNode(bfourteen,1,new TGeoCombiTrans(1488*mm/2-100.*mm/2,-100.*mm/2,0.,vinrot));
1130   vbox->AddNode(bfourteen,2,new TGeoCombiTrans(-1488*mm/2+100.*mm/2,-100.*mm/2,0.,vinrot)); 
1131   TGeoRotation *rotboxbar = new TGeoRotation("rotboxbar"); rotboxbar->RotateZ(-35);
1132   TGeoRotation *arotboxbar = new TGeoRotation("arotboxbar"); arotboxbar->RotateZ(-35); arotboxbar->RotateY(180);
1133   vbox->AddNode(dfourteen,1,new TGeoCombiTrans(-1488*mm/4,-1,0.4,rotboxbar)); 
1134   vbox->AddNode(dfourteen,2,new TGeoCombiTrans(+1488*mm/4,-1,0.4,arotboxbar));
1135  //vertical box on the short base of the cradle  
1136   vboxlast->AddNode(afourteen,1,new TGeoTranslation(0.,487.*mm/2 -100.*mm/2,0.));
1137   vboxlast->AddNode(bfourteen,1,new TGeoCombiTrans(1488*mm/2-100.*mm/2,-100.*mm/2,0.,vinrot));
1138   vboxlast->AddNode(bfourteen,2,new TGeoCombiTrans(-1488*mm/2+100.*mm/2,-100.*mm/2,0.,vinrot)); 
1139   vboxlast->AddNode(dfourteenlast,1,new TGeoCombiTrans(-1488*mm/4+1.7,-3.,0.,rotboxbar)); 
1140   vboxlast->AddNode(dfourteenlast,2,new TGeoCombiTrans(+1488*mm/4-1.7,-3.,0.,arotboxbar));
1141    
1142
1143   //POSITIONING IN THE VIRTUAL VOLUME "cradle" 
1144   
1145   //long base
1146   TGeoRotation *rotl=new TGeoRotation("Clongbase"); rotl->RotateX(90);  
1147   cradle->AddNode(lbase,0,new TGeoCombiTrans (   0*mm,   (1488-100)*mm/2, -(597-60)*mm/2,rotl)); 
1148   cradle->AddNode(lbase,1,new TGeoCombiTrans (   0*mm,   -(1488-100)*mm/2, -(597-60)*mm/2,rotl)); 
1149   //short base
1150   TGeoRotation *rots=new TGeoRotation("Cshortbase"); rots->RotateX(90); rots->RotateZ(90);
1151   cradle->AddNode(sbase,1,new TGeoCombiTrans ((6037-100)*mm/2, 0.,-(597-60)*mm/2,rots));
1152   cradle->AddNode(sbase,2,new TGeoCombiTrans (-(6037-100)*mm/2, 0.,-(597-60)*mm/2,rots));
1153   
1154   //trapezoidal structure
1155   Double_t origintrapstructure = (6037-2*60)*mm/2 - 2288*mm;
1156   
1157   TGeoRotation *rot1=new TGeoRotation("inclrot"); rot1->RotateX(90); rot1->RotateY(200);
1158   TGeoRotation *rot2=new TGeoRotation("horizrot"); rot2->RotateX(-90);
1159   Double_t dx =(1607-35)*mm*TMath::Cos(TMath::DegToRad()*20)/2-tubeheight/2*TMath::Sin(TMath::DegToRad()*20)+params[5];
1160   
1161   
1162   cradle->AddNode(inclin,1,new TGeoCombiTrans(origintrapstructure + (2288+60)*mm -dx,729*mm,params[0]+0.4,rot1));//+0.7 added
1163   cradle->AddNode(horiz,1,new TGeoCombiTrans( origintrapstructure,729*mm, 597*mm/2 - tubeheight/2,rot2));//correctly positioned
1164   TGeoRotation *rot1mirror=new TGeoRotation("inclmirrot"); rot1mirror->RotateX(90); rot1mirror->RotateY(200); rot1mirror->RotateZ(180);
1165   cradle->AddNode(inclin,2,new TGeoCombiTrans(origintrapstructure - 2345*mm + dx,729*mm,params[0]+0.4,rot1mirror));//+0.7 added
1166   cradle->AddNode(inclin,3,new TGeoCombiTrans(origintrapstructure + (2288+60)*mm -dx,-729*mm,params[0]+0.4,rot1));//0.7 added
1167   cradle->AddNode(horiz,2,new TGeoCombiTrans( origintrapstructure,-729*mm, 597*mm/2 - tubeheight/2,rot2));//correctly positioned
1168   cradle->AddNode(inclin,4,new TGeoCombiTrans(origintrapstructure - 2345*mm + dx,-729*mm,params[0]+0.4,rot1mirror));//0.7 added
1169   
1170   //inner pieces on one side
1171   TGeoRotation *rot4=new TGeoRotation("4rot"); rot4->RotateX(-90); rot4->RotateY(-55); rot4->RotateZ(180);
1172   TGeoRotation *rot4a=new TGeoRotation("4arot"); rot4a->RotateX(-90); rot4a->RotateY(-55);
1173   cradle->AddNode(four,1,new TGeoCombiTrans(origintrapstructure -(39+(597-50-60)/(2*TMath::Tan(TMath::DegToRad()*55)))*mm- tubeheight/(2*TMath::Sin(TMath::DegToRad()*55)),-729*mm,params[3],rot4));
1174   
1175   cradle->AddNode(four,2,new TGeoCombiTrans(origintrapstructure +(39+(597-50-60)/(2*TMath::Tan(TMath::DegToRad()*55)))*mm+tubeheight/(2*TMath::Sin(TMath::DegToRad()*55)),-729*mm,params[3],rot4a));
1176   
1177   TGeoRotation *rot5=new TGeoRotation("5rot"); rot5->RotateX(-90); rot5->RotateY(-75); rot5->RotateZ(180);
1178   TGeoRotation *rot5a=new TGeoRotation("5arot"); rot5a->RotateX(-90); rot5a->RotateY(-75);
1179   cradle->AddNode(five,1,new TGeoCombiTrans(origintrapstructure +(486+(597-50-60)/(2*TMath::Tan(TMath::DegToRad()*75)))*mm +tubeheight/(2*TMath::Sin(TMath::DegToRad()*75)),-729*mm,0,rot5));
1180   cradle->AddNode(five,2,new TGeoCombiTrans(origintrapstructure -(486+(597-50-60)/(2*TMath::Tan(TMath::DegToRad()*75)))*mm - tubeheight/(2*TMath::Sin(TMath::DegToRad()*75)),-729*mm,0,rot5a));
1181   cradle->AddNode(six,1,new TGeoCombiTrans(origintrapstructure+808*mm+(480*mm/2)*TMath::Cos(TMath::DegToRad()*55)+tubeheight/(2*TMath::Sin(TMath::DegToRad()*55)) +
1182   2.,-729*mm,-params[4]-0.5,rot4a));
1183   cradle->AddNode(six,2,new TGeoCombiTrans(origintrapstructure-808*mm-(480*mm/2)*TMath::Cos(TMath::DegToRad()*55)-tubeheight/(2*TMath::Sin(TMath::DegToRad()*55))
1184   -2.,-729*mm,-params[4]-0.5,rot4));
1185   
1186   TGeoRotation *rot7=new TGeoRotation("7rot"); rot7->RotateX(-90); rot7->RotateY(130); rot7->RotateZ(180);
1187   TGeoRotation *rot7a=new TGeoRotation("7arot"); rot7a->RotateX(-90); rot7a->RotateY(130);
1188   
1189   cradle->AddNode(seven,1,new TGeoCombiTrans(origintrapstructure+1478*mm-(472*mm/2)*TMath::Cos(TMath::DegToRad()*50)+tubeheight/(2*TMath::Sin(TMath::DegToRad()*50)),-729*mm,-params[8],rot7));
1190   cradle->AddNode(seven,2,new
1191                   TGeoCombiTrans(origintrapstructure-1478*mm+(472*mm/2)*TMath::Cos(TMath::DegToRad()*50)-tubeheight/(2*TMath::Sin(TMath::DegToRad()*50)),-729*mm,-params[8],rot7a));
1192   TGeoRotation *rot8=new TGeoRotation("8rot"); rot8->RotateX(-90); rot8->RotateY(-25);
1193   TGeoRotation *rot8a=new TGeoRotation("8arot"); rot8a->RotateX(-90); rot8a->RotateY(-25); rot8a->RotateZ(180);
1194   cradle->AddNode(eight,1,new TGeoCombiTrans(origintrapstructure+1640*mm+(244*mm/2)*TMath::Cos(TMath::DegToRad()*30)+tubeheight/(2*TMath::Sin(TMath::DegToRad()*30)),-729*mm,-20.5,rot8));
1195   cradle->AddNode(eight,2,new
1196                   TGeoCombiTrans(origintrapstructure-1640*mm-(244*mm/2)*TMath::Cos(TMath::DegToRad()*30)-tubeheight/(2*TMath::Sin(TMath::DegToRad()*30)),-729*mm,-20.5,rot8a));
1197   TGeoRotation *rot9=new TGeoRotation("9rot"); rot9->RotateX(-90); rot9->RotateY(-90);
1198   TGeoRotation *rot9a=new TGeoRotation("9arot"); rot9a->RotateX(-90); rot9a->RotateY(-90); rot9a->RotateZ(180);
1199   cradle->AddNode(nine,1,new TGeoCombiTrans(origintrapstructure+1960*mm+2.5+3.,-729.*mm,-20.,rot9));
1200   cradle->AddNode(nine,2,new TGeoCombiTrans(origintrapstructure-1960*mm-2.5-3.,-729.*mm,-20.,rot9a));
1201  //inner pieces on the other side  
1202   TGeoRotation *rot10=new TGeoRotation("10rot"); rot10->RotateX(-90); rot10->RotateY(-120);
1203   TGeoRotation *rot10a=new TGeoRotation("10arot"); rot10a->RotateX(-90); rot10a->RotateY(-120); rot10a->RotateZ(180);
1204
1205   cradle->AddNode(ten,1,new TGeoCombiTrans(origintrapstructure+1738*mm+tubeheight/(2*TMath::Sin(TMath::DegToRad()*60))-2,+729.*mm,-13.,rot10));
1206   cradle->AddNode(ten,2,new TGeoCombiTrans(origintrapstructure-1738*mm-tubeheight/(2*TMath::Sin(TMath::DegToRad()*60))+2,+729.*mm,-13.,rot10a));
1207  
1208   TGeoRotation *rot11=new TGeoRotation("11rot"); rot11->RotateX(-90); rot11->RotateY(50);
1209   TGeoRotation *rot11a=new TGeoRotation("11arot"); rot11a->RotateX(-90); rot11a->RotateY(50); rot11a->RotateZ(180);
1210   cradle->AddNode(eleven,1,new TGeoCombiTrans(origintrapstructure-1738*mm-tubeheight/(2*TMath::Sin(TMath::DegToRad()*60))+352.*mm,+729.*mm,-12.7,rot11));
1211   cradle->AddNode(eleven,2,new TGeoCombiTrans(origintrapstructure+1738*mm+tubeheight/(2*TMath::Sin(TMath::DegToRad()*60))-352.*mm,+729.*mm,-12.7,rot11a));
1212  
1213   TGeoRotation *rot12=new TGeoRotation("12rot"); rot12->RotateX(-90); rot12->RotateY(-120);
1214   TGeoRotation *rot12a=new TGeoRotation("12arot"); rot12a->RotateX(-90); rot12a->RotateY(-120); rot12a->RotateZ(180);
1215   cradle->AddNode(twelve,1,new TGeoCombiTrans(origintrapstructure+1065*mm,+729.*mm,1.,rot12));
1216   cradle->AddNode(twelve,2,new TGeoCombiTrans(origintrapstructure-1065*mm,+729.*mm,1.,rot12a));
1217  
1218  
1219    TGeoRotation *rot13=new TGeoRotation("13rot"); rot13->RotateX(-90); rot13->RotateY(-43); rot13->RotateZ(180);
1220    TGeoRotation *rot13a=new TGeoRotation("13arot"); rot13a->RotateX(-90); rot13a->RotateY(-43);
1221    cradle->AddNode(thirteen,1,new TGeoCombiTrans(origintrapstructure+572*mm - 18.,+729.*mm,-1.5,rot13));
1222    cradle->AddNode(thirteen,2,new TGeoCombiTrans(origintrapstructure-572*mm + 18.,+729.*mm,-1.5,rot13a));
1223
1224 //vertical structures
1225   TGeoRotation *vrot = new TGeoRotation("vertbox"); vrot->RotateX(90); vrot->RotateZ(90);
1226   cradle->AddNode(vboxlast,1,new TGeoCombiTrans(-6037*mm/2+50.*mm/2,0.,0.5,vrot));//vertial box on the short cradle base
1227   
1228   cradle->AddNode(vbox,2,new TGeoCombiTrans(-6037*mm/2+50.*mm/2+990.*mm,0.,0.5,vrot));
1229   cradle->AddNode(cfourteen,2,new TGeoCombiTrans(-6037*mm/2+50.*mm/2+990.*mm,0.,-477.*mm/2 -20.*mm/2,vrot));
1230   
1231   cradle->AddNode(vbox, 3,    new TGeoCombiTrans(origintrapstructure-(1641.36*mm+params[7])/2. + 50.*mm/2. +3, 0.,   0.5,vrot));
1232   cradle->AddNode(cfourteen,3,new TGeoCombiTrans(origintrapstructure-(1641.36*mm+params[7])/2. + 50.*mm/2. +3, 0.,-477.*mm/2 -20.*mm/2,vrot));
1233
1234   cradle->AddNode(vbox,4,new TGeoCombiTrans(origintrapstructure+(1641.36*mm+params[7])/2. - 50.*mm/2. -3,0.,0.5,vrot));
1235   cradle->AddNode(cfourteen,4,new TGeoCombiTrans(origintrapstructure+(1641.36*mm+params[7])/2. - 50.*mm/2. -3,0.,-477.*mm/2 -20.*mm/2,vrot));
1236
1237 return cradle;
1238 }//CreateCradle()
1239
1240
1241 TGeoVolume * AliHMPIDv3::CradleBaseVolume(TGeoMedium *med, Double_t l[7],const char *name)
1242 {
1243 /*
1244 The trapezoid is build in the xy plane
1245
1246     0  ________________ 1
1247       /       |        \
1248      /        |         \
1249     /       (0,0)        \
1250    /          |           \
1251 3 /___________|____________\ 2
1252
1253  01 is right shifted => shift is positive
1254
1255   //1: small base (0-1); 2: long base (3-2);
1256   //3: trapezoid height; 4: shift between the two bases;
1257   //5: height 6: height reduction; 7: z-reduction;
1258 */
1259
1260   
1261   TGeoXtru   *xtruIn  = new TGeoXtru(2);
1262   TGeoXtru   *xtruOut = new TGeoXtru(2);
1263   xtruIn->SetName(Form("%sIN",name));
1264   xtruOut->SetName(Form("%sOUT",name));
1265   
1266   Double_t xv[4], yv[4];
1267   
1268   xv[0] = -(l[0]/2 - l[3]); yv[0] =  l[2]/2;
1269   xv[1] =  l[0]/2 + l[3];   yv[1] =  l[2]/2;
1270   xv[2] =  l[1]/2; yv[2] = -l[2]/2;
1271   xv[3] = -l[1]/2; yv[3] = -l[2]/2;
1272   
1273   xtruOut->DefinePolygon(4, xv, yv);
1274   xtruOut->DefineSection(0, -l[4]/2., 0., 0., 1.0);//0=  I plane z; (0.,0.) = shift wrt centre; 1.= shape scale factor 
1275   xtruOut->DefineSection(1, +l[4]/2., 0., 0., 1.0);//1= II plane z;
1276   
1277   Double_t tgalpha = 0;
1278    if(xv[3]-xv[0] == 0 ) tgalpha = 999999; 
1279    else tgalpha =  l[2]/TMath::Abs(xv[3]-xv[0]);
1280   Double_t tgbeta = 0;
1281   if(xv[2]-xv[1]==0) tgbeta = 999999;
1282   else tgbeta = l[2]/TMath::Abs(xv[2]-xv[1]);
1283   
1284   xv[0] = xv[0]-l[5]/tgalpha; yv[0] =  l[2]/2 - l[5];
1285   xv[1] =  xv[1]+l[5]/tgbeta; yv[1] =  l[2]/2 - l[5];
1286   xv[2] =  xv[2]+l[5]/tgbeta; yv[2] = -l[2]/2+l[5];
1287   xv[3] = xv[3]-l[5]/tgalpha; yv[3] = -l[2]/2+l[5];
1288   
1289   xtruIn->DefinePolygon(4, xv, yv);
1290   xtruIn->DefineSection(0, (-l[4]+l[6])/2, 0., 0., 1.0);
1291   xtruIn->DefineSection(1,  (+l[4]-l[6])/2, 0., 0., 1.0);
1292   
1293   TGeoCompositeShape *shape = new TGeoCompositeShape(name, Form("%sOUT-%sIN",name,name));
1294   
1295   TGeoVolume *vol = new TGeoVolume(name, shape, med);
1296   
1297   return vol;
1298 }//CradleBaseVolume()