]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSvBeamTestITS04.cxx
Modifications to reduce the size in memory (E.Crescio)
[u/mrichter/AliRoot.git] / ITS / AliITSvBeamTestITS04.cxx
1 ////////////////////////////////////////////////////////
2 // ITS geometry class and step manager for the        //
3 //   integrated ITS test beam of Nov. 04              //
4 //  Author: mercedes.lopez.noriega@cern.ch            //
5 ////////////////////////////////////////////////////////
6 #include "AliRun.h"
7 #include "AliITSvBeamTestITS04.h"
8 #include <TClonesArray.h>
9 #include <TString.h>
10 #include "AliITS.h"
11 #include "AliITSDetTypeSim.h"
12 #include "AliITSgeom.h"
13 #include "AliITShit.h"
14 #include "AliITSresponseSDD.h"
15 #include "AliITSresponseSPD.h"
16 #include "AliITSresponseSSD.h"
17 #include "AliITSsegmentationSDD.h"
18 #include "AliITSsegmentationSPD.h"
19 #include "AliITSsegmentationSSD.h"
20 #include "AliMagF.h"
21 #include "TVirtualMC.h"
22 #include "AliMC.h"
23
24 const Int_t AliITSvBeamTestITS04::fgkNumberOfSPD = 4;
25 const Int_t AliITSvBeamTestITS04::fgkNumberOfSDD = 2;
26 const Int_t AliITSvBeamTestITS04::fgkNumberOfSSD = 4;
27 // Dimension (thickness:Y (beam direction), width:X, length:Z)
28
29 const char*    AliITSvBeamTestITS04::fgSPDsensitiveVolName = "ITSspdSensitiv";
30 //dimensions (preliminary values from Petra (in cms))
31 const Double_t AliITSvBeamTestITS04::fgkSPDthickness    = 0.02;
32 const Double_t AliITSvBeamTestITS04::fgkSPDwidth        = 1.4; 
33 const Double_t AliITSvBeamTestITS04::fgkSPDlength       = 7.2;
34 const Double_t AliITSvBeamTestITS04::fgkSPDthickSens    = 0.02;
35 const Double_t AliITSvBeamTestITS04::fgkSPDwidthSens    = 1.2; 
36 const Double_t AliITSvBeamTestITS04::fgkSPDlengthSens   = 7.0;
37 //position
38 const Double_t AliITSvBeamTestITS04::fgkSPD0y = 23.7;
39 const Double_t AliITSvBeamTestITS04::fgkSPD1y = 33.7;
40
41 //===
42 const char*    AliITSvBeamTestITS04::fgSDDsensitiveVolName = "ITSsddSensitiv";
43 //dimensions (preliminary values from Ludovic (in cms))
44 const Double_t AliITSvBeamTestITS04::fgkSDDthickness     = 0.03;
45 const Double_t AliITSvBeamTestITS04::fgkSDDwidth         = 7.22;
46 const Double_t AliITSvBeamTestITS04::fgkSDDlength        = 8.76;
47 const Double_t AliITSvBeamTestITS04::fgkSDDthickSens     = 0.02998;
48 const Double_t AliITSvBeamTestITS04::fgkSDDwidthSens     = 7.017;
49 const Double_t AliITSvBeamTestITS04::fgkSDDlengthSens    = 7.497;
50 //position
51 const Double_t AliITSvBeamTestITS04::fgkSDD0y = 51.7;
52 const Double_t AliITSvBeamTestITS04::fgkSDD1y = 57.2;
53
54 //===
55 const char*    AliITSvBeamTestITS04::fgSSDsensitiveVolName = "ITSssdSensitiv";
56 //dimensions (final values from Javier (in cms))
57 const Double_t AliITSvBeamTestITS04::fgkSSDthickness    = 0.03;
58 const Double_t AliITSvBeamTestITS04::fgkSSDwidth        = 7.7;
59 const Double_t AliITSvBeamTestITS04::fgkSSDlength       = 4.4;
60 const Double_t AliITSvBeamTestITS04::fgkSSDthickSens    = 0.03;
61 const Double_t AliITSvBeamTestITS04::fgkSSDwidthSens    = 7.5;
62 const Double_t AliITSvBeamTestITS04::fgkSSDlengthSens   = 4.2;
63 //position
64 const Double_t AliITSvBeamTestITS04::fgkSSD0y = 73.6;
65 const Double_t AliITSvBeamTestITS04::fgkSSD1y = 80.6;
66
67 //===============================================================
68
69
70 #include <TLorentzVector.h>
71 #include "AliTrackReference.h"
72 #include "AliITSDetTypeSim.h"
73 #include "AliITSgeom.h"
74 #include "AliITSgeomSDD.h"
75 #include "AliITSgeomSPD.h"
76 #include "AliITSgeomSSD.h"
77 #include "AliITShit.h"
78 #include "AliITSresponseSDD.h"
79 #include "AliITSresponseSPD.h"
80 #include "AliITSresponseSSD.h"
81 #include "AliITSsegmentationSDD.h"
82 #include "AliITSsegmentationSPD.h"
83 #include "AliITSsegmentationSSD.h"
84
85 #include <TGeoManager.h>
86 #include <TGeoVolume.h>
87 #include <TGeoPcon.h>
88
89 ClassImp(AliITSvBeamTestITS04)
90     
91 //_____________________________________________________________
92 AliITSvBeamTestITS04::AliITSvBeamTestITS04() : AliITS()
93 {
94     //
95     // Constructor
96     //
97     
98     SetNumberOfSPD(fgkNumberOfSPD);
99     SetNumberOfSDD(fgkNumberOfSDD);
100     SetNumberOfSSD(fgkNumberOfSSD);
101     
102     fIdN = 3;         
103     fIdName    = new TString[fIdN];
104     fIdName[0] = fgSPDsensitiveVolName;
105     fIdName[1] = fgSDDsensitiveVolName;
106     fIdName[2] = fgSSDsensitiveVolName;
107     fIdSens    = new Int_t[fIdN];
108     for(Int_t i=0; i<fIdN; i++) fIdSens[i] = 0;
109     
110     //for writing out geometry
111     fGeomDetOut   = kFALSE; 
112
113     // for reading in geometry (JC)
114     fGeomDetIn = kFALSE;
115
116     for(Int_t a=0;a<60;a++) fWrite[a] = '\0';
117 }
118
119 //_____________________________________________________________
120 AliITSvBeamTestITS04::AliITSvBeamTestITS04(const char* name,const char *title)
121     : AliITS(name,title)
122 {
123     //
124     // Constructor
125     //
126     
127     SetNumberOfSPD(fgkNumberOfSPD);
128     SetNumberOfSDD(fgkNumberOfSDD);
129     SetNumberOfSSD(fgkNumberOfSSD);
130     
131     fIdN = 3;         
132     fIdName    = new TString[fIdN];
133     fIdName[0] = fgSPDsensitiveVolName;
134     fIdName[1] = fgSDDsensitiveVolName;
135     fIdName[2] = fgSSDsensitiveVolName;
136     fIdSens    = new Int_t[fIdN];
137     for(Int_t i=0; i<fIdN; i++) fIdSens[i] = 0;
138     
139     //for writing out geometry
140     fGeomDetOut   = kFALSE; // Don't write .det file
141
142     // for reading in geometry (JC)
143     fGeomDetIn = kFALSE;
144
145     for(Int_t a=0;a<60;a++) fWrite[a] = '\0';    
146 }
147
148 //______________________________________________________________________
149 AliITSvBeamTestITS04::AliITSvBeamTestITS04(const AliITSvBeamTestITS04 &source) :  AliITS(source){
150   //Copy constructor (dummy)
151     if(&source == this) return;
152     Warning("Copy Constructor","Not allowed to copy AliITSvSDD03");
153     return;
154 }
155 //______________________________________________________________________
156 AliITSvBeamTestITS04& AliITSvBeamTestITS04::operator=(const AliITSvBeamTestITS04 &source){
157
158
159     // This class is not to be copied. Function only dummy.
160     if(&source == this) return *this;
161     Warning("= operator","Not allowed to copy AliITSvSDD03");
162     return *this;
163 }
164
165 //__________________________________________________________________
166 AliITSvBeamTestITS04::~AliITSvBeamTestITS04()
167 {
168     //
169     // Destructor
170     //
171 }
172
173 //______________________________________________________________________
174 void AliITSvBeamTestITS04::CreateMaterials()
175 {
176     // Media defined here should correspond to the one defined in galice.cuts
177     // This file is read in (AliMC*) fMCApp::Init() { ReadTransPar(); }
178     
179     // Create ITS materials
180     Int_t   ifield = gAlice->Field()->Integ();
181     Float_t fieldm = gAlice->Field()->Max();
182     
183     Float_t tmaxfdSi = 0.1;
184     Float_t stemaxSi = 0.0075;
185     Float_t deemaxSi = 0.1;
186     Float_t epsilSi  = 1.0E-4;
187     Float_t stminSi  = 0.0;
188     
189     Float_t tmaxfdAir = 0.1;
190     Float_t stemaxAir = .10000E+01;
191     Float_t deemaxAir = 0.1;
192     Float_t epsilAir  = 1.0E-4;
193     Float_t stminAir  = 0.0;
194  
195     // AIR
196     Float_t aAir[4]={12.0107,14.0067,15.9994,39.948};
197     Float_t zAir[4]={6.,7.,8.,18.};
198     Float_t wAir[4]={0.000124,0.755267,0.231781,0.012827};
199     Float_t dAir = 1.20479E-3;
200     
201     AliMaterial(51,"ITSspdSi",0.28086E+02,0.14000E+02,0.23300E+01,0.93600E+01,0.99900E+03);
202     AliMedium(51,"ITSspdSi",51,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,epsilSi,stminSi);
203     
204     AliMaterial(1,"ITSsddSi",0.28086E+02,0.14000E+02,0.23300E+01,0.93600E+01,0.99900E+03);
205     AliMedium(1,"ITSsddSi",1,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,epsilSi,stminSi);
206     
207     //AliMaterial(?,"ITSssdSi",0.28086E+02,0.14000E+02,0.23300E+01,0.93600E+01,0.99900E+03);
208     //AliMedium(?,"ITSssdSi",51,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,epsilSi,stminSi);
209     
210     AliMixture(5,"ITSair",aAir,zAir,dAir,4,wAir);
211     AliMedium(5,"ITSair",5,0,ifield,fieldm,tmaxfdAir,stemaxAir,deemaxAir,epsilAir,stminAir);
212     
213 //NEED TO ADD PLASTIC OF SCINTILLATORS!!
214
215 }
216
217 //______________________________________________________________________
218 void AliITSvBeamTestITS04::CreateGeometry()
219 {    
220   //Creates geometry
221     TGeoManager *geoManager = gGeoManager;
222     TGeoVolume *vALIC = geoManager->GetTopVolume();
223     
224     //================================
225     //       ITS mother volume
226     //================================
227     TGeoPcon *sITS = new TGeoPcon("ITS Top Volume",0.0,360.0,2);
228     // DefineSection(section number, Z, Rmin, Rmax).
229     sITS->DefineSection(0,-100.0,0.01,100.0); // Units in cms
230     sITS->DefineSection(1,+100.0,0.01,100.0);
231     
232     TGeoMedium *air = gGeoManager->GetMedium("ITSair");
233     fITSmotherVolume = new TGeoVolume("vITS",sITS,air);
234     fITSmotherVolume->SetVisibility(kFALSE);
235     vALIC->AddNode(fITSmotherVolume,1,0);
236     
237 //     //Scintillators
238 //     TGeoMedium *plasticScint = new TGeoMedium("plasticScint",1,Plastic);
239 //     //First Scintillator
240 //     TGeoBBox *Scint1Shape = new TGeoBBox("Scint1Shape",0.5,0.1,0.5,0); //1x1cm
241 //     TGeoVolume *Scint1 = new TGeoVolume("Scint1",Scint1Shape,plasticScint);
242 //     TGeoTranslation *firstScint = new TGeoTranslation(0,0.7,0);
243 //     vALIC->AddNode(Scint1,2,firstScint);
244 //     //Second Scintillator
245 //     TGeoBBox *Scint2Shape = new TGeoBBox("Scint2Shape",1.,0.1,1.,0); //2x2cm
246 //     TGeoVolume *Scint2 = new TGeoVolume("Scint2",Scint2Shape,plasticScint);
247 //     TGeoTranslation *secondScint = new TGeoTranslation(0,90.,0);
248 //     vALIC->AddNode(Scint2,3,secondScint);
249     
250     AddSPDGeometry(fITSmotherVolume);
251     AddSDDGeometry(fITSmotherVolume);
252     AddSSDGeometry(fITSmotherVolume);
253 }
254
255 //______________________________________________________________________
256 void AliITSvBeamTestITS04::Init()
257 {
258     // Initialize the ITS after it has been created.
259     Int_t i;
260     for(i=0;i<20;i++) printf("*");
261     printf( " ITSbeamtest_Init " );
262     for(i=0;i<20;i++) printf("*"); printf("\n");
263
264 //    // Create geometry
265 //    if(!fGeomDetIn) this->InitAliITSgeom();
266
267     // Initialize AliITS
268     AliITS::Init();
269     for(i=0;i<40+16;i++) printf("*"); printf("\n");
270
271 }
272
273 //______________________________________________________________________
274 void AliITSvBeamTestITS04::InitAliITSgeom()
275 {    
276   //initialisation of ITSgeom
277     const Int_t knlayers = 6;
278     Int_t nlad[knlayers], ndet[knlayers];
279     
280     nlad[0] = 1; ndet[0] = 2;
281     nlad[1] = 1; ndet[1] = 2;
282     nlad[2] = 1; ndet[2] = 1;
283     nlad[3] = 1; ndet[3] = 1;
284     nlad[4] = 1; ndet[4] = 2;
285     nlad[5] = 1; ndet[5] = 2;
286
287     Int_t nModTot = fNspd + fNsdd + fNssd;
288     if (GetITSgeom()) SetITSgeom(0x0);
289     AliITSgeom* geom = new AliITSgeom(0,knlayers,nlad,ndet,nModTot);
290     SetITSgeom(geom);
291     //*** Set default shapes 
292     const Float_t kDxyzSPD[] = {fgkSPDwidthSens/2, fgkSPDthickSens/2,fgkSPDlengthSens/2};  
293     if(!(GetITSgeom()->IsShapeDefined(kSPD)))
294         GetITSgeom()->ReSetShape(kSPD,new AliITSgeomSPD425Short(3,(Float_t *)kDxyzSPD));
295     
296     const Float_t kDxyzSDD[] = {fgkSDDwidthSens/2., fgkSDDthickSens/2.,fgkSDDlengthSens/2.};
297     if(!(GetITSgeom()->IsShapeDefined(kSDD)))
298         GetITSgeom()->ReSetShape(kSDD, new AliITSgeomSDD256(3,(Float_t *)kDxyzSDD));
299     
300     const Float_t kDxyzSSD[] = {fgkSSDlengthSens/2, fgkSSDthickSens/2,fgkSSDwidthSens/2};
301     if(!(GetITSgeom()->IsShapeDefined(kSSD)))
302         GetITSgeom()->ReSetShape(kSSD,new AliITSgeomSSD75and275(3,(Float_t *)kDxyzSSD));
303     
304     // Creating the matrices in AliITSgeom for each sensitive volume
305     // (like in AliITSv11GeometrySDD) mln
306     // Here, each layer is one detector
307     
308     char layerName[30];
309     Int_t startMod = 0;
310     
311     // SPD
312     for (Int_t i=0; i<fNspd;i++) {
313         sprintf(layerName, "ITSspdWafer_%i",i+1);
314         TGeoNode *layNode = fITSmotherVolume->GetNode(layerName);
315         if (layNode) {
316             TGeoHMatrix layMatrix(*layNode->GetMatrix());           
317             Double_t *trans  = layMatrix.GetTranslation();
318             Double_t *r      = layMatrix.GetRotationMatrix();
319             Double_t rot[10] = {r[0],r[1],r[2],
320                                 r[3],r[4],r[5],
321                                 r[6],r[7],r[8], 1.0};
322             Int_t iDet = 1;
323             if ((i+1==2)||(i+1==4)) iDet = 2;
324             Int_t iLad = 1;
325             Int_t iLay = 1;
326             if (i+1>2) iLay = 2;
327             GetITSgeom()->CreatMatrix(startMod,iLay,iLad,iDet,kSPD,trans,rot);
328             startMod++;
329         };
330     };
331     
332     // SDD
333     for (Int_t i=0; i<fNsdd;i++) {
334         sprintf(layerName, "ITSsddWafer_%i",i+fNspd+1);
335         TGeoNode *layNode = fITSmotherVolume->GetNode(layerName);
336         if (layNode) {
337             TGeoHMatrix layMatrix(*layNode->GetMatrix());
338             Double_t *trans  = layMatrix.GetTranslation();
339             Double_t *r      = layMatrix.GetRotationMatrix();
340             Double_t rot[10] = {r[0],r[1],r[2],
341                                 r[3],r[4],r[5],
342                                 r[6],r[7],r[8], 1.0};
343             Int_t iDet = 1;
344             Int_t iLad = 1;
345             Int_t iLay = fNspd-1+i;
346             GetITSgeom()->CreatMatrix(startMod,iLay,iLad,iDet,kSDD,trans,rot);
347             startMod++;
348         };
349     };
350     
351     // SSD
352     for (Int_t i=0; i<fNssd;i++) {
353         sprintf(layerName, "ITSssdWafer_%i",i+fNspd+fNsdd+1);
354         TGeoNode *layNode = fITSmotherVolume->GetNode(layerName);
355         if (layNode) {
356             TGeoHMatrix layMatrix(*layNode->GetMatrix());           
357             Double_t *trans  = layMatrix.GetTranslation();
358             Double_t *r      = layMatrix.GetRotationMatrix();
359             Double_t rot[10] = {r[0],r[1],r[2],
360                                 r[3],r[4],r[5],
361                                 r[6],r[7],r[8], 1.0};
362             Int_t iDet = 1;
363             if ((i+1==2)||(i+1==4)) iDet = 2;
364             Int_t iLad = 1;
365             Int_t iLay = 5;
366             if (i+1>2) iLay = 6;
367             GetITSgeom()->CreatMatrix(startMod,iLay,iLad,iDet,kSSD,trans,rot);
368             startMod++;
369         };
370     };
371     
372     return;
373 }
374
375 //______________________________________________________________________
376 void AliITSvBeamTestITS04::SetDefaults()
377 {
378     // (from AliITSv11) mln
379     
380     const Float_t kconv = 1.0e+04; // convert cm to microns
381     
382     if(!fDetTypeSim) fDetTypeSim = new AliITSDetTypeSim();
383     fDetTypeSim->SetITSgeom(GetITSgeom());
384     fDetTypeSim->ResetResponse();
385     fDetTypeSim->ResetSegmentation();
386  
387     AliITSgeomSPD *s0;
388     AliITSgeomSDD *s1;
389     AliITSgeomSSD *s2;
390     Int_t i;
391     Float_t bx[256],bz[280];
392
393     // If fGeomDetIn is set true the geometry will
394     // be initialised from file (JC)
395     if(GetITSgeom()!=0) SetITSgeom(0x0);
396     AliITSgeom* geom = new AliITSgeom();
397     SetITSgeom(geom);
398     if(fGeomDetIn) GetITSgeom()->ReadNewFile(fRead);
399     if(!fGeomDetIn) this->InitAliITSgeom();
400     if(fGeomDetOut) GetITSgeom()->WriteNewFile(fWrite);
401
402    
403     // SPD
404
405     s0 = (AliITSgeomSPD*) GetITSgeom()->GetShape(kSPD);// Get shape info.
406     if (s0) {
407         AliITSresponse *resp0=new AliITSresponseSPD();
408         SetResponseModel(kSPD,resp0);
409
410         AliITSsegmentationSPD *seg0=new AliITSsegmentationSPD(GetITSgeom());
411         seg0->SetDetSize(s0->GetDx()*2.*kconv, // base this on AliITSgeomSPD
412                          s0->GetDz()*2.*kconv, // for now.
413                          s0->GetDy()*2.*kconv);// x,z,y full width in microns.
414         seg0->SetNPads(256,160);               // Number of Bins in x and z
415         for(i=000;i<256;i++) bx[i] =  50.0;    // in x all are 50 microns.
416         for(i=000;i<160;i++) bz[i] = 425.0;    // most are 425 microns except below
417         for(i=160;i<280;i++) bz[i] =   0.0;    // Outside of detector.
418         bz[ 31] = bz[ 32] = 625.0;             // first chip boundry
419         bz[ 63] = bz[ 64] = 625.0;             // first chip boundry
420         bz[ 95] = bz[ 96] = 625.0;             // first chip boundry
421         bz[127] = bz[128] = 625.0;             // first chip boundry
422         bz[160] = 425.0;                       // Set so that there is no zero pixel size for fNz.
423         seg0->SetBinSize(bx,bz);               // Based on AliITSgeomSPD for now.
424         SetSegmentationModel(kSPD,seg0);
425         // set digit and raw cluster classes to be used
426         const char *kData0=(fDetTypeSim->GetResponseModel(kSPD))->DataType();
427         if (strstr(kData0,"real")) 
428           fDetTypeSim->SetDigitClassName(kSPD,"AliITSdigit");
429         else fDetTypeSim->SetDigitClassName(kSPD,"AliITSdigitSPD");
430     }
431   
432     // SDD
433    
434     s1 = (AliITSgeomSDD*) GetITSgeom()->GetShape(kSDD);// Get shape info.
435     if (s1) {
436       AliITSresponseSDD *resp1=new AliITSresponseSDD("simulated");
437       SetResponseModel(kSDD,resp1);
438       AliITSsegmentationSDD *seg1=new AliITSsegmentationSDD(GetITSgeom(),resp1);
439       seg1->SetDetSize(s1->GetDx()*kconv, // base this on AliITSgeomSDD
440                        s1->GetDz()*4.*kconv, // for now.
441                        s1->GetDy()*4.*kconv); // x,z,y full width in microns.
442       seg1->SetNPads(256,256);// Use AliITSgeomSDD for now
443       SetSegmentationModel(kSDD,seg1);
444       const char *kData1=(fDetTypeSim->GetResponseModel(kSDD))->DataType();
445       const char *kopt=fDetTypeSim->GetResponseModel(kSDD)->ZeroSuppOption();
446       if((!strstr(kopt,"2D")) && (!strstr(kopt,"1D")) || strstr(kData1,"real") ){
447         fDetTypeSim->SetDigitClassName(kSDD,"AliITSdigit");
448         } else fDetTypeSim->SetDigitClassName(kSDD,"AliITSdigitSDD");
449     }
450     
451     // SSD
452     
453     s2 = (AliITSgeomSSD*) GetITSgeom()->GetShape(kSSD);// Get shape info. Do it this way for now.
454     if (s2) {
455       AliITSresponse *resp2=new AliITSresponseSSD("simulated");
456       SetResponseModel(kSSD,resp2);
457
458       AliITSsegmentationSSD *seg2=new AliITSsegmentationSSD(GetITSgeom());
459       seg2->SetDetSize(s2->GetDx()*2.*kconv, // base this on AliITSgeomSSD
460                        s2->GetDz()*2.*kconv, // for now.
461                        s2->GetDy()*2.*kconv); // x,z,y full width in microns.
462       seg2->SetPadSize(95.,0.); // strip x pitch in microns
463       seg2->SetNPads(768,0); // number of strips on each side.
464       seg2->SetAngles(0.0075,0.0275); // strip angels rad P and N side.
465       seg2->SetAnglesLay5(0.0075,0.0275); // strip angels rad P and N side.
466       seg2->SetAnglesLay6(0.0275,0.0075); // strip angels rad P and N side.
467       SetSegmentationModel(kSSD,seg2); 
468       const char *kData2=(fDetTypeSim->GetResponseModel(kSSD))->DataType();
469       if(strstr(kData2,"real") ) fDetTypeSim->SetDigitClassName(kSSD,"AliITSdigit");
470       else fDetTypeSim->SetDigitClassName(kSSD,"AliITSdigitSSD");
471     }
472     
473   if(fgkNTYPES>3){Warning("SetDefaults","Only the four basic detector types are initialised!");}
474   return;
475 }
476
477 //______________________________________________________________________
478 void AliITSvBeamTestITS04::AddSPDGeometry(TGeoVolume *moth) const
479 {
480   //Adds SPD geometry
481     TGeoMedium *siliconSPD = gGeoManager->GetMedium("ITSspdSi");
482     
483     //outer volume
484     TGeoBBox *waferSPDshape = new TGeoBBox("ITSspdWaferShape",fgkSPDwidth/2,fgkSPDthickness/2,fgkSPDlength/2,0);
485     TGeoVolume *waferSPD = new TGeoVolume("ITSspdWafer",waferSPDshape,siliconSPD);
486     //sensitive volume
487     TGeoBBox *sensSPDbox = new TGeoBBox("ITSsddSensorSensBox",fgkSPDwidthSens/2,fgkSPDthickSens/2,fgkSPDlengthSens/2,0);
488     TGeoVolume *sensVolSPD = new TGeoVolume(fgSPDsensitiveVolName,sensSPDbox,siliconSPD);
489     waferSPD->AddNode(sensVolSPD, 1, 0); //added to outer volume
490     
491     //locate them in space (with respect top volume)
492     TGeoTranslation *spd1tr = new TGeoTranslation(0,fgkSPD0y,fgkSPDlength/2);
493     TGeoTranslation *spd2tr = new TGeoTranslation(0,fgkSPD0y,-fgkSPDlength/2);
494
495     TGeoTranslation *spd3tr = new TGeoTranslation(0,fgkSPD1y,fgkSPDlength/2);
496     TGeoTranslation *spd4tr = new TGeoTranslation(0,fgkSPD1y,-fgkSPDlength/2);
497     
498     //add to top volume
499     moth->AddNode(waferSPD, 1, spd1tr);
500     moth->AddNode(waferSPD, 2, spd2tr);
501     moth->AddNode(waferSPD, 3, spd3tr);
502     moth->AddNode(waferSPD, 4, spd4tr);
503     
504     //draw options
505     waferSPD->SetLineColor(4);
506     sensVolSPD->SetLineColor(4);
507 }
508
509
510 //______________________________________________________________________
511 void AliITSvBeamTestITS04::AddSDDGeometry(TGeoVolume *moth) const
512 {
513   //Adds SDD geometry
514     TGeoMedium *siliconSDD = gGeoManager->GetMedium("ITSsddSi");
515     
516     //outer volume
517     TGeoBBox *waferSDDshape = new TGeoBBox("ITSsddWaferShape",fgkSDDwidth/2,fgkSDDthickness/2,fgkSDDlength/2,0);
518     TGeoVolume *waferSDD = new TGeoVolume("ITSsddWafer",waferSDDshape,siliconSDD);
519     //sensitive volume
520     TGeoBBox *sensSDDbox = new TGeoBBox("ITSsddSensorSensBox",fgkSDDwidthSens/2,fgkSDDthickSens/2,fgkSDDlengthSens/2,0);
521     TGeoVolume *sensVolSDD = new TGeoVolume(fgSDDsensitiveVolName,sensSDDbox,siliconSDD);
522     waferSDD->AddNode(sensVolSDD, 1, 0); //added to outer volume
523     
524     //locate them in space
525     TGeoTranslation *sdd1tr = new TGeoTranslation(0,fgkSDD0y,0);
526     TGeoTranslation *sdd2tr = new TGeoTranslation(0,fgkSDD1y,0);
527         
528     //add to top volume
529     moth->AddNode(waferSDD, fNspd+1, sdd1tr);
530     moth->AddNode(waferSDD, fNspd+2, sdd2tr);
531     
532     //draw options
533     waferSDD->SetLineColor(3);
534     sensVolSDD->SetLineColor(3);
535 }
536
537
538 //______________________________________________________________________
539 void AliITSvBeamTestITS04::AddSSDGeometry(TGeoVolume *moth) const
540 {
541   //Adds SSD geometry
542     TGeoMedium *siliconSSD = gGeoManager->GetMedium("ITSspdSi"); // SSD medium still needed!!!
543     
544     //outer volume 
545     TGeoBBox *waferSSDshape = new TGeoBBox("ITSssdWaferShape",fgkSSDwidth/2,fgkSSDthickness/2,fgkSSDlength/2,0);
546     TGeoVolume *waferSSD = new TGeoVolume("ITSssdWafer",waferSSDshape,siliconSSD);
547     //sensitive volume
548     TGeoBBox *sensSSDbox = new TGeoBBox("ITSssdSensorSensBox",fgkSSDwidthSens/2,fgkSSDthickSens/2,fgkSSDlengthSens/2,0);
549     TGeoVolume *sensVolSSD = new TGeoVolume(fgSSDsensitiveVolName,sensSSDbox,siliconSSD);
550     waferSSD->AddNode(sensVolSSD, 1, 0);
551     
552     //locate them in space
553     /* In the SSD, there was an overlap of sensitive volumes of 2.9mm = 0.29cm (0.29/2=0.145) 
554        in the modules in the same plane, therefore the modules where not in the same plane in 
555        the Y direction, there was a "thickness" (0.03cm) difference */
556     TGeoTranslation *ssd1tr = new TGeoTranslation(0,fgkSSD0y,fgkSSDlength/2-0.145);
557     TGeoTranslation *ssd2tr = new TGeoTranslation(0,fgkSSD0y+0.03,-fgkSSDlength/2+0.145);
558
559     TGeoTranslation *ssd3tr = new TGeoTranslation(0,fgkSSD1y,fgkSSDlength/2-0.145);
560     TGeoTranslation *ssd4tr = new TGeoTranslation(0,fgkSSD1y+0.03,-fgkSSDlength/2+0.145);
561
562     //add to top volume
563     moth->AddNode(waferSSD, fNspd+fNsdd+1, ssd1tr);
564     moth->AddNode(waferSSD, fNspd+fNsdd+2, ssd2tr);
565     moth->AddNode(waferSSD, fNspd+fNsdd+3, ssd3tr);
566     moth->AddNode(waferSSD, fNspd+fNsdd+4, ssd4tr);
567     
568     //draw options
569     waferSSD->SetLineColor(2);
570     sensVolSSD->SetLineColor(2);
571 }
572
573 //______________________________________________________________________
574 void AliITSvBeamTestITS04::StepManager()
575 {
576     // Called for every step in the ITS, then calles the AliITShit class
577     // creator with the information to be recoreded about that hit.
578
579     // "Standard" StepManager. (Similar to AliITSv11) mln
580     Int_t copy, id;
581     TLorentzVector position, momentum;
582     static TLorentzVector position0;
583     static Int_t stat0=0;
584     
585     if(!(this->IsActive())){
586         return;
587     } // end if !Active volume.
588     
589     if(!(gMC->TrackCharge())) return;
590     
591     id=gMC->CurrentVolID(copy);
592     
593     Bool_t sensvol = kFALSE;
594     for(Int_t kk = 0; kk < fIdN; kk++)
595         if(id == fIdSens[kk]) sensvol = kTRUE;
596     
597     if (sensvol && (gMC->IsTrackExiting())) {
598         copy = fTrackReferences->GetEntriesFast();
599         TClonesArray &lTR = *fTrackReferences;
600         // Fill TrackReference structure with this new TrackReference.
601         new(lTR[copy]) AliTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber());
602     } // if Outer ITS mother Volume
603     
604     Int_t   vol[5];
605     TClonesArray &lhits = *fHits;
606     //
607     // Track status
608     vol[3] = 0;
609     vol[4] = 0;
610     // Fill hit structure.
611     if(gMC->IsTrackInside())      vol[3] +=  1;
612     if(gMC->IsTrackEntering())    vol[3] +=  2;
613     if(gMC->IsTrackExiting())     vol[3] +=  4;
614     if(gMC->IsTrackOut())         vol[3] +=  8;
615     if(gMC->IsTrackDisappeared()) vol[3] += 16;
616     if(gMC->IsTrackStop())        vol[3] += 32;
617     if(gMC->IsTrackAlive())       vol[3] += 64;
618     
619     // Only entering charged tracks
620     if(!(gMC->TrackCharge())) return;
621     
622     if( ((id = gMC->CurrentVolID(copy)) == fIdSens[0]) ||
623         ((id = gMC->CurrentVolID(copy)) == fIdSens[1]) ||
624         ((id = gMC->CurrentVolID(copy)) == fIdSens[2]) )
625     {
626         GetCurrentLayLaddDet(vol[0], vol[2], vol[1]);
627         // vol[2], vol[1]) : in this order because the ladder
628         // index and the det. index are exchanged in the constructor
629         // of AliITShit...
630     } else {
631         return; // not an ITS volume?
632     };
633     
634     gMC->TrackPosition(position);
635     gMC->TrackMomentum(momentum);
636     vol[4] = stat0;
637     if(gMC->IsTrackEntering()){
638         position0 = position;
639         stat0 = vol[3];
640         return;
641     } // end if IsEntering
642     // Fill hit structure with this new hit.
643     new(lhits[fNhits++]) AliITShit(fIshunt,gAlice->GetMCApp()->GetCurrentTrackNumber(),
644                                    vol, gMC->Edep(),gMC->TrackTime(),position,
645                                    position0,momentum);
646     //
647     position0 = position;
648     stat0 = vol[3];
649     return;
650 }
651
652 //______________________________________________________________________
653 Int_t AliITSvBeamTestITS04::GetCurrentLayLaddDet(Int_t &lay,Int_t &ladd, Int_t &det) const
654
655   // Function which gives the layer, ladder and det.
656   // index of the current volume. To be used in
657   // AliITS::StepManager()
658
659     det  = 1;   ladd = 1;
660     
661     TGeoNode *node = gGeoManager->GetMother(1);
662     if (!node) return kFALSE;
663     Int_t nodeNum = node->GetNumber();
664     
665     // GetNumber() return the index recorded in the node
666     
667     if (nodeNum==5||nodeNum==6) {         // SDD: one layer, one detector
668         lay = nodeNum-2;
669     } else if (nodeNum==3||nodeNum==4) {  // SPD layer 2
670         lay = 2;
671         if (nodeNum==4) det = 2;
672     } else if (nodeNum==1||nodeNum==2){   // SPD layer 1
673         lay = 1;
674         if (nodeNum==2) det = 2; 
675     } else if (nodeNum==9||nodeNum==10) { // SSD layer 2
676         lay = 6;
677         if (nodeNum==10) det = 2;
678     } else if (nodeNum==7||nodeNum==8){   // SSD layer 1
679         lay = 5;
680         if (nodeNum==8) det = 2; 
681     };  
682     
683     return kTRUE;
684 }
685
686 //_____________________________________________________________
687
688  Int_t AliITSvBeamTestITS04::GetNumberOfSubDet(const TString& det) const{
689     
690    //Get number of individual detectors
691     if(det.Contains("SPD")) return fNspd;
692     if(det.Contains("SDD")) return fNsdd;
693     if(det.Contains("SSD")) return fNssd;
694     return 0;
695   }