]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HMPID/AliHMPIDv3.cxx
Example macro for merging sets of alignment objects
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDv3.cxx
CommitLineData
100711d2 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 gMC
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()
100711d2 32#include <TF1.h> //DefineOpticalProperties()
33#include <TF2.h> //DefineOpticalProperties()
34#include <TLorentzVector.h> //IsLostByFresnel()
35#include <AliCDBManager.h> //CreateMaterials()
36#include <AliCDBEntry.h> //CreateMaterials()
37#include <TGeoPhysicalNode.h> //AddAlignableVolumes()
38#include <TGeoXtru.h> //CradleBaseVolume()
39#include <TGeoCompositeShape.h> //CradleBaseVolume()
40#include <TString.h> //StepManager()
db53fc59 41#include "AliGeomManager.h" //AddAlignableVolumes()
42
100711d2 43ClassImp(AliHMPIDv3)
44//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
45void AliHMPIDv3::AddAlignableVolumes()const
46{
47// Associates the symbolic volume name with the corresponding volume path. Interface method from AliModule invoked from AliMC
48// Arguments: none
49// Returns: none
50
db53fc59 51 AliGeomManager::ELayerID idHMPID = AliGeomManager::kHMPID;
52 Int_t modUID, modnum = 0;
53
100711d2 54 TGeoHMatrix *pGm = new TGeoHMatrix;
55 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)
56 pGm->SetTranslation(trans);
57
58 Double_t ph[7]={10.,10., 30.,30.,30. ,50.,50};
59
60 for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++) {
db53fc59 61 modUID = AliGeomManager::LayerToVolUID(idHMPID,modnum++);
62 if(!gGeoManager->SetAlignableEntry(Form("/HMPID/Chamber%i",iCh),Form("ALIC_1/Hmp%i_0",iCh),modUID))
5b89b470 63 AliError("AliHMPIDv3::Unable to set alignable entry!!"); //aligment without AliCluster3D
db53fc59 64 //Get Tracking To Local matricies for alignment with AliCluster3D
65 TGeoPNEntry *eCh = gGeoManager->GetAlignableEntryByUID(modUID);
66 TGeoHMatrix *globMatrix = eCh->GetGlobalOrig();
67
68 //Double_t phi = 20.0 * ((iCh+1) / 3) + 10.0;
69 Double_t phi = ph[iCh];
70 TGeoHMatrix *t2l = new TGeoHMatrix();
71 t2l->RotateZ(phi);
72 t2l->MultiplyLeft(&(globMatrix->Inverse()));
73 eCh->SetMatrix(t2l);
74 }//iCh loop
75
100711d2 76}
77//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
78void AliHMPIDv3::CreateMaterials()
79{
80// Definition of available HMPID materials
81// Arguments: none
82// Returns: none
83 AliDebug(1,"Start v2 HMPID.");
84
85 //clm update material definition later on from Antonello
86
87//data from PDG booklet 2002 density [gr/cm^3] rad len [cm] abs len [cm]
88 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
89 Float_t aC6F14[2]={ 12.01 , 18.99} , zC6F14[2]={ 6 , 9} , wC6F14[2]={6 , 14} , dC6F14=1.68 ; Int_t nC6F14=-2;
90 Float_t aSiO2[2]={ 28.09 , 15.99} , zSiO2[2]={14 , 8} , wSiO2[2]={1 , 2} , dSiO2=2.64 ; Int_t nSiO2=-2;
91 Float_t aCH4[2]={ 12.01 , 1.01} , zCH4[2]={ 6 , 1} , wCH4[2]={1 , 4} , dCH4=7.17e-4 ; Int_t nCH4=-2;
92// 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;
93
94 Float_t aRoha = 12.01 , zRoha = 6 , dRoha = 0.10 , radRoha = 18.80 , absRoha = 86.3/dRoha; //special material- quasi quartz
95 Float_t aCu = 63.55 , zCu = 29 , dCu = 8.96 , radCu = 1.43 , absCu = 134.9/dCu ;
96 Float_t aW =183.84 , zW = 74 , dW = 19.30 , radW = 0.35 , absW = 185.0/dW ;
97 Float_t aAl = 26.98 , zAl = 13 , dAl = 2.70 , radAl = 8.90 , absAl = 106.4/dAl ;
98 Float_t aAr = 39.94 , zAr = 18 , dAr = 1.396e-3, radAr = 14.0 , absAr = 117.2/dAr ;
99
100 Int_t matId=0; //tmp material id number
101 Int_t unsens = 0, sens=1; //sensitive or unsensitive medium
102 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
103 Float_t maxfld = gAlice->Field()->Max(); //max field value
104 Float_t tmaxfd = -10.0; //max deflection angle due to magnetic field in one step
105 Float_t deemax = - 0.2; //max fractional energy loss in one step
106 Float_t stemax = - 0.1; //max step allowed [cm]
107 Float_t epsil = 0.001; //abs tracking precision [cm]
108 Float_t stmin = - 0.001; //min step size [cm] in continius process transport, negative value: choose it automatically
109
110 // PCB copmposed mainly by G10 (Si,C,H,O) -> CsI is negligible (<500nm thick)
111 // So what is called CsI has the optical properties of CsI, but the composition of G-10 (for delta elec, etc production...)
112
113 Float_t aG10[4] = {28.09,12.01,1.01,16.00};
114 Float_t zG10[4] = {14., 6., 1., 8.};
115 Float_t wG10[4] = {0.129060,0.515016,0.061873,0.294050};
116 Float_t dG10 = 1.7;
117 Int_t nG10 = 4;
118
119 AliMixture(++matId,"Air" ,aAir ,zAir ,dAir ,nAir ,wAir ); AliMedium(kAir ,"Air" ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
120 AliMixture(++matId,"C6F14",aC6F14,zC6F14,dC6F14,nC6F14,wC6F14); AliMedium(kC6F14,"C6F14",matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
121 AliMixture(++matId,"SiO2" ,aSiO2 ,zSiO2 ,dSiO2 ,nSiO2 ,wSiO2 ); AliMedium(kSiO2 ,"SiO2" ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
122 AliMixture(++matId,"CH4" ,aCH4 ,zCH4 ,dCH4 ,nCH4 ,wCH4 ); AliMedium(kCH4 ,"CH4" ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
123// AliMixture(++matId,"CsI" ,aCsI ,zCsI ,dCsI ,nCsI ,wCsI ); AliMedium(kCsI ,"CsI" ,matId, sens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);//sensitive
124 AliMixture(++matId,"CsI+PCB",aG10 , zG10, dG10,nG10 ,wG10 ); AliMedium(kCsI ,"CsI" ,matId, sens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);//sensitive
125
126 AliMixture(++matId ,"Neo" ,aSiO2 ,zSiO2 ,dSiO2 ,nSiO2 ,wSiO2 ); AliMedium(kNeo ,"Neo" ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin); //clm neoceram
127 AliMaterial(++matId,"Roha",aRoha,zRoha,dRoha,radRoha,absRoha); AliMedium(kRoha ,"Roha" ,matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin); //Roha->honeycomb
128
129
130 AliMaterial(++matId,"Cu" ,aCu ,zCu ,dCu ,radCu ,absCu ); AliMedium(kCu ,"Cu" , matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
131 AliMaterial(++matId,"W" ,aW ,zW ,dW ,radW ,absW ); AliMedium(kW ,"W" , matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
132 AliMaterial(++matId,"Al" ,aAl ,zAl ,dAl ,radAl ,absAl ); AliMedium(kAl ,"Al" , matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
133 AliMaterial(++matId,"Ar" ,aAr ,zAr ,dAr ,radAr ,absAr ); AliMedium(kAr ,"Ar" , matId, unsens, itgfld, maxfld, tmaxfd, stemax, deemax, epsil, stmin);
134
135 InitProperties();
136
137}//void AliHMPID::CreateMaterials()
138//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
139void 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)
147HMPID 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)
149HMPID 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)
151HMPID 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)
153HMPID 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)
155HMPID 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)
157HMPID 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 gMC->Gstpar(idtmed[imed], "CUTGAM",1.e-4);
164 gMC->Gstpar(idtmed[imed], "CUTELE",1.e-4);
165 gMC->Gstpar(idtmed[imed], "CUTNEU",1.e-4);
166 gMC->Gstpar(idtmed[imed], "CUTMUO",1.e-4);
167 gMC->Gstpar(idtmed[imed], "DCUTE" ,1.e-3);
168 gMC->Gstpar(idtmed[imed], "CUTHAD",1.e-3);
169
170 gMC->Gstpar(idtmed[imed], "DRAY",1);
171 gMC->Gstpar(idtmed[imed], "LOSS",1);
172
173 imed = kC6F14; // * Freon Radiator (> 500 keV delta-electrons)
174 gMC->Gstpar(idtmed[imed], "CUTGAM",1.e-4);
175 gMC->Gstpar(idtmed[imed], "CUTELE",1.e-4);
176 gMC->Gstpar(idtmed[imed], "CUTNEU",1.e-4);
177 gMC->Gstpar(idtmed[imed], "CUTMUO",1.e-4);
178 gMC->Gstpar(idtmed[imed], "DCUTE" ,5.e-4);
179 gMC->Gstpar(idtmed[imed], "CUTHAD",5.e-4);
180
181 gMC->Gstpar(idtmed[imed], "DRAY",1);
182 gMC->Gstpar(idtmed[imed], "LOSS",1);
183
184 imed = kCH4; // * Methane Gap (> 100 keV delta-electrons)
185 gMC->Gstpar(idtmed[imed], "CUTGAM",5.e-5);
186 gMC->Gstpar(idtmed[imed], "CUTELE",5.e-5);
187 gMC->Gstpar(idtmed[imed], "CUTNEU",1.e-4);
188 gMC->Gstpar(idtmed[imed], "CUTMUO",1.e-4);
189 gMC->Gstpar(idtmed[imed], "DCUTE" ,1.e-4);
190 gMC->Gstpar(idtmed[imed], "CUTHAD",1.e-4);
191
192 gMC->Gstpar(idtmed[imed], "DRAY",1);
193 gMC->Gstpar(idtmed[imed], "LOSS",1);
194
195 imed = kCsI; // * CSI (> 50 keV delta-electrons)
196 gMC->Gstpar(idtmed[imed], "CUTGAM",1.e-5);
197 gMC->Gstpar(idtmed[imed], "CUTELE",1.e-5);
198 gMC->Gstpar(idtmed[imed], "CUTNEU",1.e-4);
199 gMC->Gstpar(idtmed[imed], "CUTMUO",1.e-4);
200 gMC->Gstpar(idtmed[imed], "DCUTE" ,5.e-5);
201 gMC->Gstpar(idtmed[imed], "CUTHAD",5.e-5);
202
203 gMC->Gstpar(idtmed[imed], "DRAY",1);
204 gMC->Gstpar(idtmed[imed], "LOSS",1);
205
206 imed = kAl; // * Alluminium (> 50 keV delta-electrons)
207 gMC->Gstpar(idtmed[imed], "CUTGAM",1.e-5);
208 gMC->Gstpar(idtmed[imed], "CUTELE",1.e-5);
209 gMC->Gstpar(idtmed[imed], "CUTNEU",1.e-4);
210 gMC->Gstpar(idtmed[imed], "CUTMUO",1.e-4);
211 gMC->Gstpar(idtmed[imed], "DCUTE" ,5.e-5);
212 gMC->Gstpar(idtmed[imed], "CUTHAD",5.e-5);
213
214 gMC->Gstpar(idtmed[imed], "DRAY",1);
215 gMC->Gstpar(idtmed[imed], "LOSS",1);
216
217 imed = kCu; // * Copper (> 50 keV delta-electrons)
218 gMC->Gstpar(idtmed[imed], "CUTGAM",1.e-5);
219 gMC->Gstpar(idtmed[imed], "CUTELE",1.e-5);
220 gMC->Gstpar(idtmed[imed], "CUTNEU",1.e-4);
221 gMC->Gstpar(idtmed[imed], "CUTMUO",1.e-4);
222 gMC->Gstpar(idtmed[imed], "DCUTE" ,5.e-5);
223 gMC->Gstpar(idtmed[imed], "CUTHAD",5.e-5);
224
225 gMC->Gstpar(idtmed[imed], "DRAY",1);
226 gMC->Gstpar(idtmed[imed], "LOSS",1);
227
228 imed = kW; // * Tungsten (> 50 keV delta-electrons)
229 gMC->Gstpar(idtmed[imed], "CUTGAM",1.e-5);
230 gMC->Gstpar(idtmed[imed], "CUTELE",1.e-5);
231 gMC->Gstpar(idtmed[imed], "CUTNEU",1.e-4);
232 gMC->Gstpar(idtmed[imed], "CUTMUO",1.e-4);
233 gMC->Gstpar(idtmed[imed], "DCUTE" ,5.e-5);
234 gMC->Gstpar(idtmed[imed], "CUTHAD",5.e-5);
235
236 gMC->Gstpar(idtmed[imed], "DRAY",1);
237 gMC->Gstpar(idtmed[imed], "LOSS",1);
238
239}
240//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
241void 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(!gMC->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
271TGeoVolume * 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
c4860469 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
100711d2 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
c4860469 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.));
100711d2 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
c4860469 416 rect->AddNode(gassipl2,1,new TGeoTranslation(0.,0.,0));
417
100711d2 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.));
100711d2 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
c4860469 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
100711d2 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
c4860469 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));
100711d2 497 fr3->AddNode( fr3up,1, new TGeoTranslation(0., 0., 7*mm));
498 fr3->AddNode(fr3down,1,new TGeoTranslation(0., 0., -10*mm));
499
500return hmp;
501
502}//CreateChamber()
503//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
504void AliHMPIDv3::Init()
505{
506// This method defines ID for sensitive volumes, i.e. such geometry volumes for which there are if(gMC->CurrentVolID()==XXX)
507// statements in StepManager()
508// Arguments: none
509// Returns: none
510 AliDebug(1,"Start v2 HMPID.");
100711d2 511 AliDebug(1,"Stop v2 HMPID.");
512}
513//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
514void AliHMPIDv3::DefineOpticalProperties()
515{
516// Optical properties definition.
517 const Int_t kNbins=30; //number of photon energy points
518 Float_t emin=5.5,emax=8.5; //Photon energy range,[eV]
519 Float_t aEckov [kNbins];
520 Double_t dEckov [kNbins];
521 Float_t aAbsRad[kNbins], aAbsWin[kNbins], aAbsGap[kNbins], aAbsMet[kNbins];
522 Float_t aIdxRad[kNbins], aIdxWin[kNbins], aIdxGap[kNbins], aIdxMet[kNbins], aIdxPc[kNbins];
523 Float_t aQeAll [kNbins], aQePc [kNbins];
524 Double_t dReflMet[kNbins], dQePc[kNbins];
525
526 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
527 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
528 TF1 *pGaIF=new TF1("HidxGap","1+0.12489e-6/(2.62e-4 - x*x/1239.84/1239.84)" ,emin,emax); //?????? from where
529
530 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
531 pRaAF->SetParameters(3.20491e16,-0.00917890,0.742402,3035.37,4.81171,0.626309);
532 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
533 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
534
535 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
536
537 TString title=GetTitle();
538 Bool_t isFlatIdx=title.Contains("FlatIdx");
539
540 for(Int_t i=0;i<kNbins;i++){
541 Float_t eV=emin+0.1*i; //Ckov energy in eV
542 aEckov [i] =1e-9*eV; //Ckov energy in GeV
543 dEckov [i] = aEckov[i];
544 aAbsRad[i]=pRaAF->Eval(eV); (isFlatIdx)? aIdxRad[i]=1.292: aIdxRad[i]=pRaIF->Eval(eV,20);
545 aAbsWin[i]=pWiAF->Eval(eV); aIdxWin[i]=pWiIF->Eval(eV);
546 aAbsGap[i]=pGaAF->Eval(eV); aIdxGap[i]=pGaIF->Eval(eV);
547 aQeAll[i] =1; //QE for all other materials except for PC must be 1.
548 aAbsMet[i] =0.0001; aIdxMet[i]=0; //metal ref idx must be 0 in order to reflect photon
549 aIdxPc [i]=1; aQePc [i]=pQeF->Eval(eV); //PC ref idx must be 1 in order to apply photon to QE conversion
550 dQePc [i]=pQeF->Eval(eV);
551 dReflMet[i] = 0.; // no reflection on the surface of the pc (?)
552 }
553 gMC->SetCerenkov((*fIdtmed)[kC6F14] , kNbins, aEckov, aAbsRad , aQeAll , aIdxRad );
554 gMC->SetCerenkov((*fIdtmed)[kSiO2] , kNbins, aEckov, aAbsWin , aQeAll , aIdxWin );
555 gMC->SetCerenkov((*fIdtmed)[kCH4] , kNbins, aEckov, aAbsGap , aQeAll , aIdxGap );
556 gMC->SetCerenkov((*fIdtmed)[kCu] , kNbins, aEckov, aAbsMet , aQeAll , aIdxMet );
557 gMC->SetCerenkov((*fIdtmed)[kW] , kNbins, aEckov, aAbsMet , aQeAll , aIdxMet ); //n=0 means reflect photons
558 gMC->SetCerenkov((*fIdtmed)[kCsI] , kNbins, aEckov, aAbsMet , aQePc , aIdxPc ); //n=1 means convert photons
559 gMC->SetCerenkov((*fIdtmed)[kAl] , kNbins, aEckov, aAbsMet , aQeAll , aIdxMet );
560
561 // Define a skin surface for the photocatode to enable 'detection' in G4
562 gMC->DefineOpSurface("surfPc", kGlisur /*kUnified*/,kDielectric_metal,kPolished, 0.);
563 gMC->SetMaterialProperty("surfPc", "EFFICIENCY", kNbins, dEckov, dQePc);
564 gMC->SetMaterialProperty("surfPc", "REFLECTIVITY", kNbins, dEckov, dReflMet);
565 gMC->SetSkinSurface("skinPc", "Rpc", "surfPc");
566
567 delete pRaAF;delete pWiAF;delete pGaAF; delete pRaIF; delete pWiIF; delete pGaIF; delete pQeF;
568}
569//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
570Bool_t AliHMPIDv3::IsLostByFresnel()
571{
572// Calculate probability for the photon to be lost by Fresnel reflection.
573 TLorentzVector p4;
574 Double_t mom[3],localMom[3];
575 gMC->TrackMomentum(p4); mom[0]=p4(1); mom[1]=p4(2); mom[2]=p4(3);
576 localMom[0]=0; localMom[1]=0; localMom[2]=0;
577 gMC->Gmtod(mom,localMom,2);
578 Double_t localTc = localMom[0]*localMom[0]+localMom[2]*localMom[2];
579 Double_t localTheta = TMath::ATan2(TMath::Sqrt(localTc),localMom[1]);
580 Double_t cotheta = TMath::Abs(TMath::Cos(localTheta));
581 if(gMC->GetRandom()->Rndm() < Fresnel(p4.E()*1e9,cotheta,1)){
582 AliDebug(1,"Photon lost");
583 return kTRUE;
584 }else
585 return kFALSE;
586}//IsLostByFresnel()
587//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
588void AliHMPIDv3::GenFee(Float_t qtot)
589{
590// Generate FeedBack photons for the current particle. To be invoked from StepManager().
591// eloss=0 means photon so only pulse height distribution is to be analysed.
592 TLorentzVector x4;
593 gMC->TrackPosition(x4);
594 Int_t iNphotons=gMC->GetRandom()->Poisson(0.02*qtot); //# of feedback photons is proportional to the charge of hit
595 AliDebug(1,Form("N photons=%i",iNphotons));
596 Int_t j;
597 Float_t cthf, phif, enfp = 0, sthf, e1[3], e2[3], e3[3], vmod, uswop,dir[3], phi,pol[3], mom[4];
598//Generate photons
599 for(Int_t i=0;i<iNphotons;i++){//feedbacks loop
600 Double_t ranf[2];
601 gMC->GetRandom()->RndmArray(2,ranf); //Sample direction
602 cthf=ranf[0]*2-1.0;
603 if(cthf<0) continue;
604 sthf = TMath::Sqrt((1 - cthf) * (1 + cthf));
605 phif = ranf[1] * 2 * TMath::Pi();
606
607 if(Double_t randomNumber=gMC->GetRandom()->Rndm()<=0.57)
608 enfp = 7.5e-9;
609 else if(randomNumber<=0.7)
610 enfp = 6.4e-9;
611 else
612 enfp = 7.9e-9;
613
614
615 dir[0] = sthf * TMath::Sin(phif); dir[1] = cthf; dir[2] = sthf * TMath::Cos(phif);
616 gMC->Gdtom(dir, mom, 2);
617 mom[0]*=enfp; mom[1]*=enfp; mom[2]*=enfp;
618 mom[3] = TMath::Sqrt(mom[0]*mom[0]+mom[1]*mom[1]+mom[2]*mom[2]);
619
620 // Polarisation
621 e1[0]= 0; e1[1]=-dir[2]; e1[2]= dir[1];
622 e2[0]=-dir[1]; e2[1]= dir[0]; e2[2]= 0;
623 e3[0]= dir[1]; e3[1]= 0; e3[2]=-dir[0];
624
625 vmod=0;
626 for(j=0;j<3;j++) vmod+=e1[j]*e1[j];
627 if (!vmod) for(j=0;j<3;j++) {
628 uswop=e1[j];
629 e1[j]=e3[j];
630 e3[j]=uswop;
631 }
632 vmod=0;
633 for(j=0;j<3;j++) vmod+=e2[j]*e2[j];
634 if (!vmod) for(j=0;j<3;j++) {
635 uswop=e2[j];
636 e2[j]=e3[j];
637 e3[j]=uswop;
638 }
639
640 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;
641 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;
642
643 phi = gMC->GetRandom()->Rndm()* 2 * TMath::Pi();
644 for(j=0;j<3;j++) pol[j]=e1[j]*TMath::Sin(phi)+e2[j]*TMath::Cos(phi);
645 gMC->Gdtom(pol, pol, 2);
646 Int_t outputNtracksStored;
647 gAlice->GetMCApp()->PushTrack(1, //transport
648 gAlice->GetMCApp()->GetCurrentTrackNumber(),//parent track
649 50000051, //PID
650 mom[0],mom[1],mom[2],mom[3], //track momentum
651 x4.X(),x4.Y(),x4.Z(),x4.T(), //track origin
652 pol[0],pol[1],pol[2], //polarization
653 kPFeedBackPhoton, //process ID
654 outputNtracksStored, //on return how many new photons stored on stack
655 1.0); //weight
656 }//feedbacks loop
657 AliDebug(1,"Stop.");
658}//GenerateFeedbacks()
659//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
660void AliHMPIDv3::Hits2SDigits()
661{
662// Interface method ivoked from AliSimulation to create a list of sdigits corresponding to list of hits. Every hit generates one or more sdigits.
663// Arguments: none
664// Returns: none
665 AliDebug(1,"Start.");
666 for(Int_t iEvt=0;iEvt < GetLoader()->GetRunLoader()->GetNumberOfEvents();iEvt++){ //events loop
667 GetLoader()->GetRunLoader()->GetEvent(iEvt); //get next event
668
669 if(!GetLoader()->TreeH()) {GetLoader()->LoadHits(); }
670 if(!GetLoader()->TreeS()) {GetLoader()->MakeTree("S"); MakeBranch("S");}//to
671
672 for(Int_t iEnt=0;iEnt<GetLoader()->TreeH()->GetEntries();iEnt++){//prims loop
673 GetLoader()->TreeH()->GetEntry(iEnt);
674 Hit2Sdi(Hits(),SdiLst());
675 }//prims loop
676 GetLoader()->TreeS()->Fill();
677 GetLoader()->WriteSDigits("OVERWRITE");
678 SdiReset();
679 }//events loop
680 GetLoader()->UnloadHits();
681 GetLoader()->UnloadSDigits();
682 AliDebug(1,"Stop.");
683}//Hits2SDigits()
684//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
685void AliHMPIDv3::Hit2Sdi(TClonesArray *pHitLst,TClonesArray *pSdiLst)
686{
687// Converts list of hits to list of sdigits.
688// Arguments: pHitLst - list of hits provided not empty
689// pSDigLst - list of sdigits where to store the results
690// Returns: none
691 for(Int_t iHit=0;iHit<pHitLst->GetEntries();iHit++){ //hits loop
692 AliHMPIDHit *pHit=(AliHMPIDHit*)pHitLst->At(iHit); //get pointer to current hit
693 pHit->Hit2Sdi(pSdiLst); //convert this hit to list of sdigits
694 }//hits loop loop
695}//Hits2Sdi()
696//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
697void AliHMPIDv3::Digits2Raw()
698{
699// Interface method invoked by AliSimulation to create raw data streams from digits. Events loop is done in AliSimulation
700// Arguments: none
701// Returns: none
702 AliDebug(1,"Start.");
703 GetLoader()->LoadDigits();
704 TTree * treeD = GetLoader()->TreeD();
705 if(!treeD) {
706 AliError("No digits tree!");
707 return;
708 }
709 treeD->GetEntry(0);
710
711
712 AliHMPIDRawStream *pRS=0x0;
713 pRS->WriteRaw(DigLst());
714
715 GetLoader()->UnloadDigits();
716 AliDebug(1,"Stop.");
717}//Digits2Raw()
718//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
719Float_t AliHMPIDv3::Fresnel(Float_t ene,Float_t pdoti, Bool_t pola)
720{
721// Correction for Fresnel ???????????
722// Arguments: ene - photon energy [GeV],
723// PDOTI=COS(INC.ANG.), PDOTR=COS(POL.PLANE ROT.ANG.)
724// Returns:
725 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,
726 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,
727 7.8,7.9,8.0,8.1,8.2,8.3,8.4,8.5};
728 Float_t csin[36] = {2.14,2.21,2.33,2.48,2.76,2.97,2.99,2.59,2.81,3.05,
729 2.86,2.53,2.55,2.66,2.79,2.96,3.18,3.05,2.84,2.81,2.38,2.11,
730 2.01,2.13,2.39,2.73,3.08,3.15,2.95,2.73,2.56,2.41,2.12,1.95,
731 1.72,1.53};
732 Float_t csik[36] = {0.,0.,0.,0.,0.,0.196,0.408,0.208,0.118,0.49,0.784,0.543,
733 0.424,0.404,0.371,0.514,0.922,1.102,1.139,1.376,1.461,1.253,0.878,
734 0.69,0.612,0.649,0.824,1.347,1.571,1.678,1.763,1.857,1.824,1.824,
735 1.714,1.498};
736 Float_t xe=ene;
737 Int_t j=Int_t(xe*10)-49;
738 Float_t cn=csin[j]+((csin[j+1]-csin[j])/0.1)*(xe-en[j]);
739 Float_t ck=csik[j]+((csik[j+1]-csik[j])/0.1)*(xe-en[j]);
740
741 //FORMULAE FROM HANDBOOK OF OPTICS, 33.23 OR
742 //W.R. HUNTER, J.O.S.A. 54 (1964),15 , J.O.S.A. 55(1965),1197
743
744 Float_t sinin=TMath::Sqrt(1-pdoti*pdoti);
745 Float_t tanin=sinin/pdoti;
746
747 Float_t c1=cn*cn-ck*ck-sinin*sinin;
748 Float_t c2=4*cn*cn*ck*ck;
749 Float_t aO=TMath::Sqrt(0.5*(TMath::Sqrt(c1*c1+c2)+c1));
750 Float_t b2=0.5*(TMath::Sqrt(c1*c1+c2)-c1);
751
752 Float_t rs=((aO-pdoti)*(aO-pdoti)+b2)/((aO+pdoti)*(aO+pdoti)+b2);
753 Float_t rp=rs*((aO-sinin*tanin)*(aO-sinin*tanin)+b2)/((aO+sinin*tanin)*(aO+sinin*tanin)+b2);
754
755
756 //CORRECTION FACTOR FOR SURFACE ROUGHNESS
757 //B.J. STAGG APPLIED OPTICS, 30(1991),4113
758
759 Float_t sigraf=18.;
760 Float_t lamb=1240/ene;
761 Float_t fresn;
762
763 Float_t rO=TMath::Exp(-(4*TMath::Pi()*pdoti*sigraf/lamb)*(4*TMath::Pi()*pdoti*sigraf/lamb));
764
765 if(pola)
766 {
767 Float_t pdotr=0.8; //DEGREE OF POLARIZATION : 1->P , -1->S
768 fresn=0.5*(rp*(1+pdotr)+rs*(1-pdotr));
769 }
770 else
771 fresn=0.5*(rp+rs);
772
773 fresn = fresn*rO;
774 return fresn;
775}//Fresnel()
776//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
777void AliHMPIDv3::Print(Option_t *option)const
778{
779// Debug printout
780 TObject::Print(option);
781}//void AliHMPID::Print(Option_t *option)const
782//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
783
784Bool_t AliHMPIDv3::Raw2SDigits(AliRawReader *pRR)
785{
786// Arguments: pRR- raw reader
787// Returns: kTRUE on success (currently ignored in AliSimulation::ConvertRaw2SDigits())
788//AliHMPIDDigit sdi; //tmp sdigit, raw digit will be converted to it
789
790 if(!GetLoader()->TreeS()) {MakeTree("S"); MakeBranch("S");}
791
792 TClonesArray *pSdiLst=SdiLst(); Int_t iSdiCnt=0; //tmp list of sdigits for all chambers
793 AliHMPIDRawStream stream(pRR);
794 while(stream.Next())
795 {
796 for(Int_t iPad=0;iPad<stream.GetNPads();iPad++) {
797 AliHMPIDDigit sdi(stream.GetPadArray()[iPad],stream.GetChargeArray()[iPad]);
798 new((*pSdiLst)[iSdiCnt++]) AliHMPIDDigit(sdi); //add this digit to the tmp list
799 }
800 }
801
802 GetLoader()->TreeS()->Fill(); GetLoader()->WriteSDigits("OVERWRITE");//write out sdigits
803 SdiReset();
804 return kTRUE;
805
806}//Raw2SDigits
807
808//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
809void AliHMPIDv3::StepCount()
810{
811// Count number of ckovs created
812}
813//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
814void AliHMPIDv3::StepHistory()
815{
816// This methode is invoked from StepManager() in order to print out
817 static Int_t iStepN;
818 const char *sParticle;
819 switch(gMC->TrackPid()){
820 case kProton: sParticle="PROTON" ;break;
821 case kNeutron: sParticle="neutron" ;break;
822 case kGamma: sParticle="gamma" ;break;
823 case 50000050: sParticle="CKOV" ;break;
824 case kPi0: sParticle="Pi0" ;break;
825 case kPiPlus: sParticle="Pi+" ;break;
826 case kPiMinus: sParticle="Pi-" ;break;
827 case kElectron: sParticle="electron" ;break;
828 default: sParticle="not known" ;break;
829 }
830
831 TString flag="fanny combination";
832 if(gMC->IsTrackAlive())
833 if(gMC->IsTrackEntering()) flag="enters to";
834 else if(gMC->IsTrackExiting()) flag="exits from";
835 else if(gMC->IsTrackInside()) flag="inside";
836 else
837 if(gMC->IsTrackStop()) flag="stoped in";
838
839 Int_t vid=0,copy=0;
840 TString path=gMC->CurrentVolName(); path.Prepend("-");path.Prepend(gMC->CurrentVolOffName(1));//current volume and his mother are always there
841 vid=gMC->CurrentVolOffID(2,copy); if(vid) {path.Prepend("-");path.Prepend(gMC->VolName(vid));}
842 vid=gMC->CurrentVolOffID(3,copy); if(vid) {path.Prepend("-");path.Prepend(gMC->VolName(vid));}
843
844
845 Printf("Step %i: %s (%i) %s %s m=%.6f GeV q=%.1f dEdX=%.4f Etot=%.4f",iStepN,sParticle,gMC->TrackPid(),flag.Data(),path.Data(),gMC->TrackMass(),gMC->TrackCharge(),gMC->Edep()*1e9,gMC->Etot());
846
847 Double_t gMcTrackPos[3]; gMC->TrackPosition(gMcTrackPos[0],gMcTrackPos[1],gMcTrackPos[2]);
848 Double_t gMcTrackPosLoc[3]; gMC->Gmtod(gMcTrackPos,gMcTrackPosLoc,1);
849 Printf("gMC 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]);
850
851
852
853 Printf("Step %i: tid=%i flags alive=%i disap=%i enter=%i exit=%i inside=%i out=%i stop=%i new=%i",
854 iStepN, gAlice->GetMCApp()->GetCurrentTrackNumber(),
855 gMC->IsTrackAlive(), gMC->IsTrackDisappeared(),gMC->IsTrackEntering(), gMC->IsTrackExiting(),
856 gMC->IsTrackInside(),gMC->IsTrackOut(), gMC->IsTrackStop(), gMC->IsNewTrack());
857
858 Float_t a,z,den,rad,abs; a=z=den=rad=abs=-1;
859 Int_t mid=gMC->CurrentMaterial(a,z,den,rad,abs);
860 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);
861
862 TArrayI proc; gMC->StepProcesses(proc);
863 Printf("Processes in this step:");
864 for ( int i = 0 ; i < proc.GetSize(); i++)
865 {
866 Printf("%s",TMCProcessName[proc.At(i)]);
867 }
868 Printf("End process list");
869
870 iStepN++;
871}//StepHistory()
872//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
873void AliHMPIDv3::StepManager()
874{
875// Full Step Manager.
876// Arguments: none
877// Returns: none
878// StepHistory(); return; //uncomment to print tracks history
879 // StepCount(); return; //uncomment to count photons
880
881 TString volname = gMC->CurrentVolName();
882
883//Treat photons
884 if((gMC->TrackPid()==50000050||gMC->TrackPid()==50000051)&&volname.Contains("Hpad")){ //photon (Ckov or feedback) hits on module PC (Hpad)
885 if(gMC->Edep()>0){ //photon survided QE test i.e. produces electron
886 if(IsLostByFresnel()){ gMC->StopTrack(); return;} //photon lost due to fersnel reflection on PC
887 Int_t tid= gMC->GetStack()->GetCurrentTrackNumber(); //take TID
888 Int_t pid= gMC->TrackPid(); //take PID
889 Float_t etot= gMC->Etot(); //total hpoton energy, [GeV]
890 Double_t x[3]; gMC->TrackPosition(x[0],x[1],x[2]); //take MARS position at entrance to PC
891 TString tmpname = volname; tmpname.Remove(0,4); Int_t idch = tmpname.Atoi(); //retrieve the chamber number
892 Float_t xl,yl; AliHMPIDParam::Instance()->Mars2Lors(idch,x,xl,yl); //take LORS position
893 new((*fHits)[fNhits++])AliHMPIDHit(idch,etot,pid,tid,xl,yl,x); //HIT for photon, position at P, etot will be set to Q
894 if(fDoFeed) GenFee(etot); //generate feedback photons etot is modified in hit ctor to Q of hit
895 }//photon hit PC and DE >0
896 }//photon hit PC
897
898
899 //Treat charged particles
900 static Float_t eloss; //need to store mip parameters between different steps
901 static Double_t in[3];
902
903 if(gMC->IsTrackEntering() && gMC->TrackCharge() && volname.Contains("Hpad")) //Trackref stored when entering in the pad volume
904 AddTrackReference(gMC->GetStack()->GetCurrentTrackNumber(), AliTrackReference::kHMPID); //for acceptance calculations
905
906
907 if(gMC->TrackCharge() && volname.Contains("Hcel")){ //charged particle in amplification gap (Hcel)
908 if(gMC->IsTrackEntering()||gMC->IsNewTrack()) { //entering or newly created
909 eloss=0; //reset Eloss collector
910 gMC->TrackPosition(in[0],in[1],in[2]); //take position at the entrance
911 }else if(gMC->IsTrackExiting()||gMC->IsTrackStop()||gMC->IsTrackDisappeared()){ //exiting or disappeared
912 eloss +=gMC->Edep(); //take into account last step Eloss
913 Int_t tid= gMC->GetStack()->GetCurrentTrackNumber(); //take TID
914 Int_t pid= gMC->TrackPid(); //take PID
915 Double_t out[3]; gMC->TrackPosition(out[0],out[1],out[2]); //take MARS position at exit
916 out[0]=0.5*(out[0]+in[0]); //
917 out[1]=0.5*(out[1]+in[1]); //take hit position at the anod plane
918 out[2]=0.5*(out[2]+in[2]);
919 TString tmpname = volname; tmpname.Remove(0,4); Int_t idch = tmpname.Atoi(); //retrieve the chamber number
920 Float_t xl,yl;AliHMPIDParam::Instance()->Mars2Lors(idch,out,xl,yl); //take LORS position
3e3bfa9b 921 if(eloss>0) {
922 new((*fHits)[fNhits++])AliHMPIDHit(idch,eloss,pid,tid,xl,yl,out); //HIT for MIP, position near anod plane, eloss will be set to Q
923 if(fDoFeed) GenFee(eloss); //generate feedback photons
924 }
100711d2 925 }else //just going inside
926 eloss += gMC->Edep(); //collect this step eloss
927 }//MIP in GAP
928
929}//StepManager()
930//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
931void AliHMPIDv3::TestPoint(Int_t ch,Float_t x,Float_t y)
932{
933// Utility method to check the validity of geometry by poviding some crucial points
934// Arguments: ch,x,y- crucial point definition (cm) in LORS
935// Returns: none
936 Double_t mars[3];
937 AliHMPIDParam::Instance()->Lors2Mars(ch,x,y,mars);
938 Printf("(ch=%i,locX=%.2f,locY=%.2f) %s",ch,x,y,gGeoManager->FindNode(mars[0],mars[1],mars[2])->GetName());
939}//TestPoint()
940//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
941void AliHMPIDv3::TestGeom()
942{
943//
944// Test method to check geometry
945//
946 //TGeoManager::Import("misaligned_geometry.root");
947 TGeoManager::Import("geometry.root");
948 for(Int_t ch=AliHMPIDParam::kMinCh;ch<=AliHMPIDParam::kMaxCh;ch++)
949 TestPoint(ch,0,0);
950}//TestPoint()
951//+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
952
953void AliHMPIDv3::IdealPosition(Int_t iCh,TGeoHMatrix *pMatrix) //ideal position of given chamber
954{
955// Construct ideal position matrix for a given chamber
956// Arguments: iCh- chamber ID; pMatrix- pointer to precreated unity matrix where to store the results
957// Returns: none
958 const Double_t kAngHor=19.5; // horizontal angle between chambers 19.5 grad
959 const Double_t kAngVer=20; // vertical angle between chambers 20 grad
960 const Double_t kAngCom=30; // common HMPID rotation with respect to x axis 30 grad
961 const Double_t kTrans[3]={490,0,0}; // center of the chamber is on window-gap surface
962 pMatrix->RotateY(90); // rotate around y since initial position is in XY plane -> now in YZ plane
963 pMatrix->SetTranslation(kTrans); // now plane in YZ is shifted along x
964 switch(iCh){
965 case 0: pMatrix->RotateY(kAngHor); pMatrix->RotateZ(-kAngVer); break; //right and down
966 case 1: pMatrix->RotateZ(-kAngVer); break; //down
967 case 2: pMatrix->RotateY(kAngHor); break; //right
968 case 3: break; //no rotation
969 case 4: pMatrix->RotateY(-kAngHor); break; //left
970 case 5: pMatrix->RotateZ(kAngVer); break; //up
971 case 6: pMatrix->RotateY(-kAngHor); pMatrix->RotateZ(kAngVer); break; //left and up
972 }
973 pMatrix->RotateZ(kAngCom); //apply common rotation in XY plane
974}
975
976
977void AliHMPIDv3::IdealPositionCradle(Int_t iCh,TGeoHMatrix *pMatrix) //ideal position of given one module of the cradle
978{
979// Construct ideal position matrix for a given module cradle
980// Arguments: iCh- chamber ID; pMatrix- pointer to precreated unity matrix where to store the results
981// Returns: none
982 const Double_t kAngHor=19.5; // horizontal angle between chambers 19.5 grad
983 const Double_t kAngVer=20; // vertical angle between chambers 20 grad
984 const Double_t kAngCom=30; // common HMPID rotation with respect to x axis 30 grad
985 const Double_t kTrans[3]={423.+ 29,0,67}; // z-center of the cradle module
986 pMatrix->RotateY(90); // rotate around y since initial position is in XY plane -> now in YZ plane
987 pMatrix->SetTranslation(kTrans); // now plane in YZ is shifted along x
988 switch(iCh){
989 case 0: pMatrix->RotateY(kAngHor); pMatrix->RotateZ(-kAngVer); break; //right and down
990 case 1: pMatrix->RotateZ(-kAngVer); break; //down
991 case 2: pMatrix->RotateY(kAngHor); break; //right
992 case 3: break; //no rotation
993 case 4: pMatrix->RotateY(-kAngHor); break; //left
994 case 5: pMatrix->RotateZ(kAngVer); break; //up
995 case 6: pMatrix->RotateY(-kAngHor); pMatrix->RotateZ(kAngVer); break; //left and up
996 }
997 pMatrix->RotateZ(kAngCom); //apply common rotation in XY plane
998}
999
1000
1001
1002TGeoVolume* AliHMPIDv3::CreateCradle()
1003{
1004
1005//Method that builds the Cradle geometry
1006//according to the base topology created
1007//in CradleBaseVolume(...)
1008
1009 Double_t mm = 0.1;
1010
1011 Double_t params[10]={0.5,10.,24.,-1,5.2,1.5,3.5,8.5,3.8,0.};
1012 TGeoMedium *med =gGeoManager->GetMedium("HMPID_Al");
1013 TGeoVolume *cradle=new TGeoVolumeAssembly("Hcradle");
1014
c4860469 1015 //Double_t baselong[7]={6037*mm-2*60*mm, 6037*mm-2*60*mm,60*mm,0.,100*mm,10*mm,10*mm};//2CRE2112P3
1016 Double_t baselong[7]={6037*mm-2*100*mm, 6037*mm-2*100*mm,60*mm,0.,100*mm,10*mm,10*mm};//2CRE2112P3
100711d2 1017 TGeoVolume *lbase = CradleBaseVolume(med,baselong,"cradleLbase");
1018 lbase->SetLineColor(kGray);
1019
1020 Double_t baseshort[7]={1288.*mm+2*100*mm, 1288.*mm+2*100*mm,60*mm,0.,100*mm,10*mm,10*mm};//2CRE2112P3
1021 TGeoVolume *sbase = CradleBaseVolume(med,baseshort,"cradleSbase");
1022 sbase->SetLineColor(kGray);
1023
1024 //one side
1025
1026 Double_t height = 30.*mm; //30 = 2*(1488/2-729) (2CRE2112P3)
1027 Double_t tubeheight = 50.*mm; Double_t heightred = 5.*mm; Double_t zred = 5.*mm;
1028 Double_t oneshift = tubeheight/TMath::Tan(TMath::DegToRad()*20.)+(1458.-35)*mm/2 - (1607-35)*mm/2;
c4860469 1029 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)
100711d2 1030 TGeoVolume *inclin = CradleBaseVolume(med,linclined,"inclinedbar");
1031 inclin->SetLineColor(kGray);
1032 Double_t lhorizontal[7] = {1641.36*mm+params[7],1659.*mm+params[7],tubeheight,0, height ,heightred,zred};
1033 TGeoVolume *horiz = CradleBaseVolume(med,lhorizontal,"horizontalbar");
1034 horiz->SetLineColor(kGray);
1035
1036 //inner bars, they are named as the numbering in 2CRE2112P3
c4860469 1037 Double_t fourshift = tubeheight/TMath::Tan(TMath::DegToRad()*55.);
1038 Double_t lfour[7] = {592*mm,592*mm,tubeheight,fourshift,height,heightred,zred};
100711d2 1039 TGeoVolume *four = CradleBaseVolume(med,lfour,"bar4");
1040 four->SetLineColor(kGray);
1041
1042 Double_t fiveshift = tubeheight/TMath::Tan(TMath::DegToRad()*75);
c4860469 1043 Double_t lfive[7] = {500.*mm,500.*mm,tubeheight,fiveshift,height,heightred,zred};
100711d2 1044 TGeoVolume *five = CradleBaseVolume(med,lfive,"bar5");
1045 five->SetLineColor(kGray);
1046
c4860469 1047 Double_t sixshift = tubeheight/TMath::Tan(TMath::DegToRad()*55)+459*mm/2-480*mm/2;
1048 Double_t lsix[7] = {456*mm,477*mm,tubeheight,sixshift,height,heightred,zred};
100711d2 1049 TGeoVolume *six = CradleBaseVolume(med,lsix,"bar6");
1050 six->SetLineColor(kGray);
1051
1052 Double_t sevenshift = tubeheight/TMath::Tan(TMath::DegToRad()*50)+472*mm/2-429.*mm/2;
1053 Double_t lseven[7] = {429*mm,472*mm,tubeheight,sevenshift,height,heightred,zred};
1054 TGeoVolume *seven = CradleBaseVolume(med,lseven,"bar7");
1055 seven->SetLineColor(kGray);
1056
1057 Double_t eightshift = tubeheight/TMath::Tan(TMath::DegToRad()*30)+244.*mm/2-200.*mm/2 -3;
1058 Double_t leight[7] = {200.*mm,244.*mm,tubeheight,eightshift,height,heightred,zred};
1059 TGeoVolume *eight = CradleBaseVolume(med,leight,"bar8");
1060 eight->SetLineColor(kGray);
1061
1062 Double_t nineshift = -tubeheight/TMath::Tan(TMath::DegToRad()*71)+83.*mm/2-66.*mm/2;
c4860469 1063 Double_t lnine[7] = {59.5*mm,76.5*mm,tubeheight,nineshift,height,heightred,zred};
100711d2 1064 TGeoVolume *nine = CradleBaseVolume(med,lnine,"bar9");
1065 nine->SetLineColor(kGray);
1066
1067 Double_t tenshift = (-tubeheight/TMath::Tan(TMath::DegToRad()*60) -221.*mm/2+195.*mm/2);
1068 Double_t lten[7] = {195.*mm,221.*mm,tubeheight,tenshift,height,heightred,zred};
1069 TGeoVolume *ten = CradleBaseVolume(med,lten,"bar10");
1070 ten->SetLineColor(kGray);
1071
1072 Double_t elevenshift = (-tubeheight/TMath::Tan(TMath::DegToRad()*70) -338.*mm/2+315.*mm/2);
c4860469 1073 Double_t leleven[7] = {308.*mm,331.*mm,tubeheight,elevenshift,height,heightred,zred};
100711d2 1074 TGeoVolume *eleven = CradleBaseVolume(med,leleven,"bar11");
1075 eleven->SetLineColor(kGray);
1076
1077 Double_t twelveshift = (-tubeheight/TMath::Tan(TMath::DegToRad()*60) -538.*mm/2+508.*mm/2);
c4860469 1078 Double_t ltwelve[7] = {507.*mm,537.*mm,tubeheight,twelveshift,height,heightred,zred};
100711d2 1079 TGeoVolume *twelve = CradleBaseVolume(med,ltwelve,"bar12");
1080 twelve->SetLineColor(kGray);
1081
c4860469 1082 Double_t thirteenshift = tubeheight/TMath::Tan(TMath::DegToRad()*43);
1083 Double_t lthirteen[7] = {708.*mm,708.*mm,tubeheight,thirteenshift,height,heightred,zred};
100711d2 1084 TGeoVolume *thirteen = CradleBaseVolume(med,lthirteen,"bar13");
1085 thirteen->SetLineColor(kGray);
1086
1087
1088 //vertical rectangles
c4860469 1089 TGeoVolume *vbox= new TGeoVolumeAssembly("Hvbox");
100711d2 1090 vbox->SetLineColor(kViolet);
1091 Double_t width = 50.*mm;
c4860469 1092
1093 TGeoVolume *vboxlast= new TGeoVolumeAssembly("Hvboxlast");//vertical structure on the short base
1094 vboxlast->SetLineColor(kViolet);
100711d2 1095
1096 Double_t barheight = 100.*mm;
1097 Double_t lAfourteen[7] = {1488.*mm,1488.*mm,barheight,0,width,heightred,zred};
1098 TGeoVolume *afourteen = CradleBaseVolume(med,lAfourteen,"bar14top");
1099 afourteen->SetLineColor(kGray);
1100
c4860469 1101 Double_t lBfourteen[7] = {387*mm,387.*mm,barheight,0,width,heightred,zred};
100711d2 1102 TGeoVolume *bfourteen = CradleBaseVolume(med,lBfourteen,"bar14vert");
1103 bfourteen->SetLineColor(kGray);
1104
1105 Double_t lCfourteen[7] = {1288.*mm,1288.*mm,barheight,0,width,heightred,zred};
1106 TGeoVolume *cfourteen = CradleBaseVolume(med,lCfourteen,"bar14bot");
1107 cfourteen->SetLineColor(kGray);
1108
c4860469 1109 Double_t oblshift = 50.*mm/ TMath::Tan(TMath::DegToRad()*35);
1110 Double_t lDfourteen[7] = {603.*mm,603.*mm,50.*mm,oblshift,width,heightred,zred};
100711d2 1111 TGeoVolume *dfourteen = CradleBaseVolume(med,lDfourteen,"bar14incl");
1112 dfourteen->SetLineColor(kGray);
1113
c4860469 1114
1115 Double_t lDfourteenlast[7] = {667.*mm,667.*mm,50.*mm,oblshift,width,heightred,zred};
1116 TGeoVolume *dfourteenlast = CradleBaseVolume(med,lDfourteenlast,"bar14incllast");
1117 dfourteenlast->SetLineColor(kGray);
1118
100711d2 1119 vbox->AddNode(afourteen,1,new TGeoTranslation(0.,487.*mm/2 -100.*mm/2,0.));
1120 TGeoRotation *vinrot = new TGeoRotation("vertbar"); vinrot->RotateZ(90);
1121 vbox->AddNode(bfourteen,1,new TGeoCombiTrans(1488*mm/2-100.*mm/2,-100.*mm/2,0.,vinrot));
1122 vbox->AddNode(bfourteen,2,new TGeoCombiTrans(-1488*mm/2+100.*mm/2,-100.*mm/2,0.,vinrot));
1123 TGeoRotation *rotboxbar = new TGeoRotation("rotboxbar"); rotboxbar->RotateZ(-35);
1124 TGeoRotation *arotboxbar = new TGeoRotation("arotboxbar"); arotboxbar->RotateZ(-35); arotboxbar->RotateY(180);
c4860469 1125 vbox->AddNode(dfourteen,1,new TGeoCombiTrans(-1488*mm/4,-1,0.4,rotboxbar));
1126 vbox->AddNode(dfourteen,2,new TGeoCombiTrans(+1488*mm/4,-1,0.4,arotboxbar));
1127 //vertical box on the short base of the cradle
1128 vboxlast->AddNode(afourteen,1,new TGeoTranslation(0.,487.*mm/2 -100.*mm/2,0.));
1129 vboxlast->AddNode(bfourteen,1,new TGeoCombiTrans(1488*mm/2-100.*mm/2,-100.*mm/2,0.,vinrot));
1130 vboxlast->AddNode(bfourteen,2,new TGeoCombiTrans(-1488*mm/2+100.*mm/2,-100.*mm/2,0.,vinrot));
1131 vboxlast->AddNode(dfourteenlast,1,new TGeoCombiTrans(-1488*mm/4+1.7,-3.,0.,rotboxbar));
1132 vboxlast->AddNode(dfourteenlast,2,new TGeoCombiTrans(+1488*mm/4-1.7,-3.,0.,arotboxbar));
100711d2 1133
1134
1135 //POSITIONING IN THE VIRTUAL VOLUME "cradle"
1136
1137 //long base
1138 TGeoRotation *rotl=new TGeoRotation("Clongbase"); rotl->RotateX(90);
c4860469 1139 cradle->AddNode(lbase,0,new TGeoCombiTrans ( 0*mm, (1488-100)*mm/2, -(597-60)*mm/2,rotl));
1140 cradle->AddNode(lbase,1,new TGeoCombiTrans ( 0*mm, -(1488-100)*mm/2, -(597-60)*mm/2,rotl));
100711d2 1141 //short base
1142 TGeoRotation *rots=new TGeoRotation("Cshortbase"); rots->RotateX(90); rots->RotateZ(90);
1143 cradle->AddNode(sbase,1,new TGeoCombiTrans ((6037-100)*mm/2, 0.,-(597-60)*mm/2,rots));
1144 cradle->AddNode(sbase,2,new TGeoCombiTrans (-(6037-100)*mm/2, 0.,-(597-60)*mm/2,rots));
1145
1146 //trapezoidal structure
1147 Double_t origintrapstructure = (6037-2*60)*mm/2 - 2288*mm;
1148
1149 TGeoRotation *rot1=new TGeoRotation("inclrot"); rot1->RotateX(90); rot1->RotateY(200);
1150 TGeoRotation *rot2=new TGeoRotation("horizrot"); rot2->RotateX(-90);
1151 Double_t dx =(1607-35)*mm*TMath::Cos(TMath::DegToRad()*20)/2-tubeheight/2*TMath::Sin(TMath::DegToRad()*20)+params[5];
1152
1153
c4860469 1154 cradle->AddNode(inclin,1,new TGeoCombiTrans(origintrapstructure + (2288+60)*mm -dx,729*mm,params[0]+0.4,rot1));//+0.7 added
100711d2 1155 cradle->AddNode(horiz,1,new TGeoCombiTrans( origintrapstructure,729*mm, 597*mm/2 - tubeheight/2,rot2));//correctly positioned
1156 TGeoRotation *rot1mirror=new TGeoRotation("inclmirrot"); rot1mirror->RotateX(90); rot1mirror->RotateY(200); rot1mirror->RotateZ(180);
c4860469 1157 cradle->AddNode(inclin,2,new TGeoCombiTrans(origintrapstructure - 2345*mm + dx,729*mm,params[0]+0.4,rot1mirror));//+0.7 added
1158 cradle->AddNode(inclin,3,new TGeoCombiTrans(origintrapstructure + (2288+60)*mm -dx,-729*mm,params[0]+0.4,rot1));//0.7 added
100711d2 1159 cradle->AddNode(horiz,2,new TGeoCombiTrans( origintrapstructure,-729*mm, 597*mm/2 - tubeheight/2,rot2));//correctly positioned
c4860469 1160 cradle->AddNode(inclin,4,new TGeoCombiTrans(origintrapstructure - 2345*mm + dx,-729*mm,params[0]+0.4,rot1mirror));//0.7 added
100711d2 1161
1162 //inner pieces on one side
1163 TGeoRotation *rot4=new TGeoRotation("4rot"); rot4->RotateX(-90); rot4->RotateY(-55); rot4->RotateZ(180);
1164 TGeoRotation *rot4a=new TGeoRotation("4arot"); rot4a->RotateX(-90); rot4a->RotateY(-55);
1165 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));
1166
1167 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));
1168
1169 TGeoRotation *rot5=new TGeoRotation("5rot"); rot5->RotateX(-90); rot5->RotateY(-75); rot5->RotateZ(180);
1170 TGeoRotation *rot5a=new TGeoRotation("5arot"); rot5a->RotateX(-90); rot5a->RotateY(-75);
1171 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));
1172 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));
1173 cradle->AddNode(six,1,new TGeoCombiTrans(origintrapstructure+808*mm+(480*mm/2)*TMath::Cos(TMath::DegToRad()*55)+tubeheight/(2*TMath::Sin(TMath::DegToRad()*55)) +
1174 2.,-729*mm,-params[4]-0.5,rot4a));
1175 cradle->AddNode(six,2,new TGeoCombiTrans(origintrapstructure-808*mm-(480*mm/2)*TMath::Cos(TMath::DegToRad()*55)-tubeheight/(2*TMath::Sin(TMath::DegToRad()*55))
1176 -2.,-729*mm,-params[4]-0.5,rot4));
1177
1178 TGeoRotation *rot7=new TGeoRotation("7rot"); rot7->RotateX(-90); rot7->RotateY(130); rot7->RotateZ(180);
1179 TGeoRotation *rot7a=new TGeoRotation("7arot"); rot7a->RotateX(-90); rot7a->RotateY(130);
1180
1181 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));
1182 cradle->AddNode(seven,2,new
1183 TGeoCombiTrans(origintrapstructure-1478*mm+(472*mm/2)*TMath::Cos(TMath::DegToRad()*50)-tubeheight/(2*TMath::Sin(TMath::DegToRad()*50)),-729*mm,-params[8],rot7a));
1184 TGeoRotation *rot8=new TGeoRotation("8rot"); rot8->RotateX(-90); rot8->RotateY(-25);
1185 TGeoRotation *rot8a=new TGeoRotation("8arot"); rot8a->RotateX(-90); rot8a->RotateY(-25); rot8a->RotateZ(180);
1186 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));
1187 cradle->AddNode(eight,2,new
1188 TGeoCombiTrans(origintrapstructure-1640*mm-(244*mm/2)*TMath::Cos(TMath::DegToRad()*30)-tubeheight/(2*TMath::Sin(TMath::DegToRad()*30)),-729*mm,-20.5,rot8a));
1189 TGeoRotation *rot9=new TGeoRotation("9rot"); rot9->RotateX(-90); rot9->RotateY(-90);
1190 TGeoRotation *rot9a=new TGeoRotation("9arot"); rot9a->RotateX(-90); rot9a->RotateY(-90); rot9a->RotateZ(180);
1191 cradle->AddNode(nine,1,new TGeoCombiTrans(origintrapstructure+1960*mm+2.5+3.,-729.*mm,-20.,rot9));
1192 cradle->AddNode(nine,2,new TGeoCombiTrans(origintrapstructure-1960*mm-2.5-3.,-729.*mm,-20.,rot9a));
100711d2 1193 //inner pieces on the other side
1194 TGeoRotation *rot10=new TGeoRotation("10rot"); rot10->RotateX(-90); rot10->RotateY(-120);
1195 TGeoRotation *rot10a=new TGeoRotation("10arot"); rot10a->RotateX(-90); rot10a->RotateY(-120); rot10a->RotateZ(180);
1196
c4860469 1197 cradle->AddNode(ten,1,new TGeoCombiTrans(origintrapstructure+1738*mm+tubeheight/(2*TMath::Sin(TMath::DegToRad()*60))-2,+729.*mm,-13.,rot10));
1198 cradle->AddNode(ten,2,new TGeoCombiTrans(origintrapstructure-1738*mm-tubeheight/(2*TMath::Sin(TMath::DegToRad()*60))+2,+729.*mm,-13.,rot10a));
100711d2 1199
1200 TGeoRotation *rot11=new TGeoRotation("11rot"); rot11->RotateX(-90); rot11->RotateY(50);
1201 TGeoRotation *rot11a=new TGeoRotation("11arot"); rot11a->RotateX(-90); rot11a->RotateY(50); rot11a->RotateZ(180);
c4860469 1202 cradle->AddNode(eleven,1,new TGeoCombiTrans(origintrapstructure-1738*mm-tubeheight/(2*TMath::Sin(TMath::DegToRad()*60))+352.*mm,+729.*mm,-12.7,rot11));
1203 cradle->AddNode(eleven,2,new TGeoCombiTrans(origintrapstructure+1738*mm+tubeheight/(2*TMath::Sin(TMath::DegToRad()*60))-352.*mm,+729.*mm,-12.7,rot11a));
100711d2 1204
1205 TGeoRotation *rot12=new TGeoRotation("12rot"); rot12->RotateX(-90); rot12->RotateY(-120);
1206 TGeoRotation *rot12a=new TGeoRotation("12arot"); rot12a->RotateX(-90); rot12a->RotateY(-120); rot12a->RotateZ(180);
1207 cradle->AddNode(twelve,1,new TGeoCombiTrans(origintrapstructure+1065*mm,+729.*mm,1.,rot12));
1208 cradle->AddNode(twelve,2,new TGeoCombiTrans(origintrapstructure-1065*mm,+729.*mm,1.,rot12a));
1209
1210
1211 TGeoRotation *rot13=new TGeoRotation("13rot"); rot13->RotateX(-90); rot13->RotateY(-43); rot13->RotateZ(180);
1212 TGeoRotation *rot13a=new TGeoRotation("13arot"); rot13a->RotateX(-90); rot13a->RotateY(-43);
1213 cradle->AddNode(thirteen,1,new TGeoCombiTrans(origintrapstructure+572*mm - 18.,+729.*mm,-1.5,rot13));
1214 cradle->AddNode(thirteen,2,new TGeoCombiTrans(origintrapstructure-572*mm + 18.,+729.*mm,-1.5,rot13a));
1215
1216//vertical structures
1217 TGeoRotation *vrot = new TGeoRotation("vertbox"); vrot->RotateX(90); vrot->RotateZ(90);
c4860469 1218 cradle->AddNode(vboxlast,1,new TGeoCombiTrans(-6037*mm/2+50.*mm/2,0.,0.5,vrot));//vertial box on the short cradle base
100711d2 1219
1220 cradle->AddNode(vbox,2,new TGeoCombiTrans(-6037*mm/2+50.*mm/2+990.*mm,0.,0.5,vrot));
1221 cradle->AddNode(cfourteen,2,new TGeoCombiTrans(-6037*mm/2+50.*mm/2+990.*mm,0.,-477.*mm/2 -20.*mm/2,vrot));
1222
1223 cradle->AddNode(vbox, 3, new TGeoCombiTrans(origintrapstructure-(1641.36*mm+params[7])/2. + 50.*mm/2. +3, 0., 0.5,vrot));
1224 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));
1225
1226 cradle->AddNode(vbox,4,new TGeoCombiTrans(origintrapstructure+(1641.36*mm+params[7])/2. - 50.*mm/2. -3,0.,0.5,vrot));
1227 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));
1228
1229return cradle;
1230}//CreateCradle()
1231
1232
1233TGeoVolume * AliHMPIDv3::CradleBaseVolume(TGeoMedium *med, Double_t l[7],char *name)
1234{
1235/*
1236The trapezoid is build in the xy plane
1237
1238 0 ________________ 1
1239 / | \
1240 / | \
1241 / (0,0) \
1242 / | \
12433 /___________|____________\ 2
1244
1245 01 is right shifted => shift is positive
1246
1247 //1: small base (0-1); 2: long base (3-2);
1248 //3: trapezoid height; 4: shift between the two bases;
1249 //5: height 6: height reduction; 7: z-reduction;
1250*/
1251
1252
1253 TGeoXtru *xtruIn = new TGeoXtru(2);
1254 TGeoXtru *xtruOut = new TGeoXtru(2);
1255 xtruIn->SetName(Form("%sIN",name));
1256 xtruOut->SetName(Form("%sOUT",name));
1257
1258 Double_t xv[4], yv[4];
1259
1260 xv[0] = -(l[0]/2 - l[3]); yv[0] = l[2]/2;
1261 xv[1] = l[0]/2 + l[3]; yv[1] = l[2]/2;
1262 xv[2] = l[1]/2; yv[2] = -l[2]/2;
1263 xv[3] = -l[1]/2; yv[3] = -l[2]/2;
100711d2 1264
1265 xtruOut->DefinePolygon(4, xv, yv);
d38002cf 1266 xtruOut->DefineSection(0, -l[4]/2., 0., 0., 1.0);//0= I plane z; (0.,0.) = shift wrt centre; 1.= shape scale factor
1267 xtruOut->DefineSection(1, +l[4]/2., 0., 0., 1.0);//1= II plane z;
100711d2 1268
1269 Double_t tgalpha = 0;
9c5fe8da 1270 if(xv[3]-xv[0] == 0 ) tgalpha = 999999;
100711d2 1271 else tgalpha = l[2]/TMath::Abs(xv[3]-xv[0]);
1272 Double_t tgbeta = 0;
9c5fe8da 1273 if(xv[2]-xv[1]==0) tgbeta = 999999;
100711d2 1274 else tgbeta = l[2]/TMath::Abs(xv[2]-xv[1]);
1275
1276 xv[0] = xv[0]-l[5]/tgalpha; yv[0] = l[2]/2 - l[5];
1277 xv[1] = xv[1]+l[5]/tgbeta; yv[1] = l[2]/2 - l[5];
1278 xv[2] = xv[2]+l[5]/tgbeta; yv[2] = -l[2]/2+l[5];
1279 xv[3] = xv[3]-l[5]/tgalpha; yv[3] = -l[2]/2+l[5];
1280
1281 xtruIn->DefinePolygon(4, xv, yv);
1282 xtruIn->DefineSection(0, (-l[4]+l[6])/2, 0., 0., 1.0);
1283 xtruIn->DefineSection(1, (+l[4]-l[6])/2, 0., 0., 1.0);
1284
1285 TGeoCompositeShape *shape = new TGeoCompositeShape(name, Form("%sOUT-%sIN",name,name));
1286
1287 TGeoVolume *vol = new TGeoVolume(name, shape, med);
1288
1289 return vol;
1290}//CradleBaseVolume()