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