]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSvBeamTestITS04.cxx
merging RecPoints and ClustersV2. All ClusterFinders produce AliITSRecPoints objects...
[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 "AliITSCalibrationSDD.h"
16 #include "AliITSCalibrationSPD.h"
17 #include "AliITSCalibrationSSD.h"
18 #include "AliITSsegmentationSDD.h"
19 #include "AliITSsegmentationSPD.h"
20 #include "AliITSsegmentationSSD.h"
21 #include "AliMagF.h"
22 #include "TVirtualMC.h"
23 #include "AliMC.h"
24
25 const Int_t AliITSvBeamTestITS04::fgkNumberOfSPD = 4;
26 const Int_t AliITSvBeamTestITS04::fgkNumberOfSDD = 2;
27 const Int_t AliITSvBeamTestITS04::fgkNumberOfSSD = 4;
28 // Dimension (thickness:Y (beam direction), width:X, length:Z)
29
30 const char*    AliITSvBeamTestITS04::fgSPDsensitiveVolName = "ITSspdSensitiv";
31 //dimensions (preliminary values from Petra (in cms))
32 const Double_t AliITSvBeamTestITS04::fgkSPDthickness    = 0.02;
33 const Double_t AliITSvBeamTestITS04::fgkSPDwidth        = 1.4; 
34 const Double_t AliITSvBeamTestITS04::fgkSPDlength       = 7.2;
35 const Double_t AliITSvBeamTestITS04::fgkSPDthickSens    = 0.02;
36 const Double_t AliITSvBeamTestITS04::fgkSPDwidthSens    = 1.2; 
37 const Double_t AliITSvBeamTestITS04::fgkSPDlengthSens   = 7.0;
38 //position
39 const Double_t AliITSvBeamTestITS04::fgkSPD0y = 23.7;
40 const Double_t AliITSvBeamTestITS04::fgkSPD1y = 33.7;
41
42 //===
43 const char*    AliITSvBeamTestITS04::fgSDDsensitiveVolName = "ITSsddSensitiv";
44 //dimensions (preliminary values from Ludovic (in cms))
45 const Double_t AliITSvBeamTestITS04::fgkSDDthickness     = 0.03;
46 const Double_t AliITSvBeamTestITS04::fgkSDDwidth         = 7.22;
47 const Double_t AliITSvBeamTestITS04::fgkSDDlength        = 8.76;
48 const Double_t AliITSvBeamTestITS04::fgkSDDthickSens     = 0.02998;
49 const Double_t AliITSvBeamTestITS04::fgkSDDwidthSens     = 7.017;
50 const Double_t AliITSvBeamTestITS04::fgkSDDlengthSens    = 7.497;
51 //position
52 const Double_t AliITSvBeamTestITS04::fgkSDD0y = 51.7;
53 const Double_t AliITSvBeamTestITS04::fgkSDD1y = 57.2;
54
55 //===
56 const char*    AliITSvBeamTestITS04::fgSSDsensitiveVolName = "ITSssdSensitiv";
57 //dimensions (final values from Javier (in cms))
58 const Double_t AliITSvBeamTestITS04::fgkSSDthickness    = 0.03;
59 const Double_t AliITSvBeamTestITS04::fgkSSDwidth        = 7.7;
60 const Double_t AliITSvBeamTestITS04::fgkSSDlength       = 4.4;
61 const Double_t AliITSvBeamTestITS04::fgkSSDthickSens    = 0.03;
62 const Double_t AliITSvBeamTestITS04::fgkSSDwidthSens    = 7.5;
63 const Double_t AliITSvBeamTestITS04::fgkSSDlengthSens   = 4.2;
64 //position
65 const Double_t AliITSvBeamTestITS04::fgkSSD0y = 73.6;
66 const Double_t AliITSvBeamTestITS04::fgkSSD1y = 80.6;
67
68 //===============================================================
69
70
71 #include <TLorentzVector.h>
72 #include "AliTrackReference.h"
73 #include "AliITSDetTypeSim.h"
74 #include "AliITSgeom.h"
75 #include "AliITSgeomSDD.h"
76 #include "AliITSgeomSPD.h"
77 #include "AliITSgeomSSD.h"
78 #include "AliITShit.h"
79 #include "AliITSCalibrationSDD.h"
80 #include "AliITSCalibrationSPD.h"
81 #include "AliITSCalibrationSSD.h"
82 #include "AliITSsegmentationSDD.h"
83 #include "AliITSsegmentationSPD.h"
84 #include "AliITSsegmentationSSD.h"
85
86 #include <TGeoManager.h>
87 #include <TGeoVolume.h>
88 #include <TGeoPcon.h>
89
90 ClassImp(AliITSvBeamTestITS04)
91     
92 //_____________________________________________________________
93 AliITSvBeamTestITS04::AliITSvBeamTestITS04() : AliITS()
94 {
95     //
96     // Constructor
97     //
98     
99     SetNumberOfSPD(fgkNumberOfSPD);
100     SetNumberOfSDD(fgkNumberOfSDD);
101     SetNumberOfSSD(fgkNumberOfSSD);
102     
103     fIdN = 3;         
104     fIdName    = new TString[fIdN];
105     fIdName[0] = fgSPDsensitiveVolName;
106     fIdName[1] = fgSDDsensitiveVolName;
107     fIdName[2] = fgSSDsensitiveVolName;
108     fIdSens    = new Int_t[fIdN];
109     for(Int_t i=0; i<fIdN; i++) fIdSens[i] = 0;
110     
111     //for writing out geometry
112     fGeomDetOut   = kFALSE; 
113
114     // for reading in geometry (JC)
115     fGeomDetIn = kFALSE;
116
117     for(Int_t a=0;a<60;a++) fWrite[a] = '\0';
118 }
119
120 //_____________________________________________________________
121 AliITSvBeamTestITS04::AliITSvBeamTestITS04(const char* name,const char *title)
122     : AliITS(name,title)
123 {
124     //
125     // Constructor
126     //
127     
128     SetNumberOfSPD(fgkNumberOfSPD);
129     SetNumberOfSDD(fgkNumberOfSDD);
130     SetNumberOfSSD(fgkNumberOfSSD);
131     
132     fIdN = 3;         
133     fIdName    = new TString[fIdN];
134     fIdName[0] = fgSPDsensitiveVolName;
135     fIdName[1] = fgSDDsensitiveVolName;
136     fIdName[2] = fgSSDsensitiveVolName;
137     fIdSens    = new Int_t[fIdN];
138     for(Int_t i=0; i<fIdN; i++) fIdSens[i] = 0;
139     
140     //for writing out geometry
141     fGeomDetOut   = kFALSE; // Don't write .det file
142
143     // for reading in geometry (JC)
144     fGeomDetIn = kFALSE;
145
146     for(Int_t a=0;a<60;a++) fWrite[a] = '\0';    
147 }
148
149 //______________________________________________________________________
150 AliITSvBeamTestITS04::AliITSvBeamTestITS04(const AliITSvBeamTestITS04 &source) :  AliITS(source){
151   //Copy constructor (dummy)
152     if(&source == this) return;
153     Warning("Copy Constructor","Not allowed to copy AliITSvSDD03");
154     return;
155 }
156 //______________________________________________________________________
157 AliITSvBeamTestITS04& AliITSvBeamTestITS04::operator=(const AliITSvBeamTestITS04 &source){
158
159
160     // This class is not to be copied. Function only dummy.
161     if(&source == this) return *this;
162     Warning("= operator","Not allowed to copy AliITSvSDD03");
163     return *this;
164 }
165
166 //__________________________________________________________________
167 AliITSvBeamTestITS04::~AliITSvBeamTestITS04()
168 {
169     //
170     // Destructor
171     //
172 }
173
174 //______________________________________________________________________
175 void AliITSvBeamTestITS04::CreateMaterials()
176 {
177     // Media defined here should correspond to the one defined in galice.cuts
178     // This file is read in (AliMC*) fMCApp::Init() { ReadTransPar(); }
179     
180     // Create ITS materials
181     Int_t   ifield = gAlice->Field()->Integ();
182     Float_t fieldm = gAlice->Field()->Max();
183     
184     Float_t tmaxfdSi = 0.1;
185     Float_t stemaxSi = 0.0075;
186     Float_t deemaxSi = 0.1;
187     Float_t epsilSi  = 1.0E-4;
188     Float_t stminSi  = 0.0;
189     
190     Float_t tmaxfdAir = 0.1;
191     Float_t stemaxAir = .10000E+01;
192     Float_t deemaxAir = 0.1;
193     Float_t epsilAir  = 1.0E-4;
194     Float_t stminAir  = 0.0;
195  
196     // AIR
197     Float_t aAir[4]={12.0107,14.0067,15.9994,39.948};
198     Float_t zAir[4]={6.,7.,8.,18.};
199     Float_t wAir[4]={0.000124,0.755267,0.231781,0.012827};
200     Float_t dAir = 1.20479E-3;
201     
202     AliMaterial(51,"ITSspdSi",0.28086E+02,0.14000E+02,0.23300E+01,0.93600E+01,0.99900E+03);
203     AliMedium(51,"ITSspdSi",51,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,epsilSi,stminSi);
204     
205     AliMaterial(1,"ITSsddSi",0.28086E+02,0.14000E+02,0.23300E+01,0.93600E+01,0.99900E+03);
206     AliMedium(1,"ITSsddSi",1,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,epsilSi,stminSi);
207     
208     //AliMaterial(?,"ITSssdSi",0.28086E+02,0.14000E+02,0.23300E+01,0.93600E+01,0.99900E+03);
209     //AliMedium(?,"ITSssdSi",51,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,epsilSi,stminSi);
210     
211     AliMixture(5,"ITSair",aAir,zAir,dAir,4,wAir);
212     AliMedium(5,"ITSair",5,0,ifield,fieldm,tmaxfdAir,stemaxAir,deemaxAir,epsilAir,stminAir);
213     
214 //NEED TO ADD PLASTIC OF SCINTILLATORS!!
215
216 }
217
218 //______________________________________________________________________
219 void AliITSvBeamTestITS04::CreateGeometry()
220 {    
221   //Creates geometry
222     TGeoManager *geoManager = gGeoManager;
223     TGeoVolume *vALIC = geoManager->GetTopVolume();
224     
225     //================================
226     //       ITS mother volume
227     //================================
228     TGeoPcon *sITS = new TGeoPcon("ITS Top Volume",0.0,360.0,2);
229     // DefineSection(section number, Z, Rmin, Rmax).
230     sITS->DefineSection(0,-100.0,0.01,100.0); // Units in cms
231     sITS->DefineSection(1,+100.0,0.01,100.0);
232     
233     TGeoMedium *air = gGeoManager->GetMedium("ITSair");
234     fITSmotherVolume = new TGeoVolume("vITS",sITS,air);
235     fITSmotherVolume->SetVisibility(kFALSE);
236     vALIC->AddNode(fITSmotherVolume,1,0);
237     
238 //     //Scintillators
239 //     TGeoMedium *plasticScint = new TGeoMedium("plasticScint",1,Plastic);
240 //     //First Scintillator
241 //     TGeoBBox *Scint1Shape = new TGeoBBox("Scint1Shape",0.5,0.1,0.5,0); //1x1cm
242 //     TGeoVolume *Scint1 = new TGeoVolume("Scint1",Scint1Shape,plasticScint);
243 //     TGeoTranslation *firstScint = new TGeoTranslation(0,0.7,0);
244 //     vALIC->AddNode(Scint1,2,firstScint);
245 //     //Second Scintillator
246 //     TGeoBBox *Scint2Shape = new TGeoBBox("Scint2Shape",1.,0.1,1.,0); //2x2cm
247 //     TGeoVolume *Scint2 = new TGeoVolume("Scint2",Scint2Shape,plasticScint);
248 //     TGeoTranslation *secondScint = new TGeoTranslation(0,90.,0);
249 //     vALIC->AddNode(Scint2,3,secondScint);
250     
251     AddSPDGeometry(fITSmotherVolume);
252     AddSDDGeometry(fITSmotherVolume);
253     AddSSDGeometry(fITSmotherVolume);
254 }
255
256 //______________________________________________________________________
257 void AliITSvBeamTestITS04::Init()
258 {
259     // Initialize the ITS after it has been created.
260     Int_t i;
261     for(i=0;i<20;i++) printf("*");
262     printf( " ITSbeamtest_Init " );
263     for(i=0;i<20;i++) printf("*"); printf("\n");
264
265 //    // Create geometry
266 //    if(!fGeomDetIn) this->InitAliITSgeom();
267
268     // Initialize AliITS
269     AliITS::Init();
270     for(i=0;i<40+16;i++) printf("*"); printf("\n");
271
272 }
273
274 //______________________________________________________________________
275 void AliITSvBeamTestITS04::InitAliITSgeom()
276 {    
277   //initialisation of ITSgeom
278     const Int_t knlayers = 6;
279     Int_t nlad[knlayers], ndet[knlayers];
280     
281     nlad[0] = 1; ndet[0] = 2;
282     nlad[1] = 1; ndet[1] = 2;
283     nlad[2] = 1; ndet[2] = 1;
284     nlad[3] = 1; ndet[3] = 1;
285     nlad[4] = 1; ndet[4] = 2;
286     nlad[5] = 1; ndet[5] = 2;
287
288     Int_t nModTot = fNspd + fNsdd + fNssd;
289     if (GetITSgeom()) SetITSgeom(0x0);
290     AliITSgeom* geom = new AliITSgeom(0,knlayers,nlad,ndet,nModTot);
291     SetITSgeom(geom);
292     //*** Set default shapes 
293     const Float_t kDxyzSPD[] = {fgkSPDwidthSens/2, fgkSPDthickSens/2,fgkSPDlengthSens/2};  
294     if(!(GetITSgeom()->IsShapeDefined(kSPD)))
295         GetITSgeom()->ReSetShape(kSPD,new AliITSgeomSPD425Short(3,(Float_t *)kDxyzSPD));
296     
297     const Float_t kDxyzSDD[] = {fgkSDDwidthSens/2., fgkSDDthickSens/2.,fgkSDDlengthSens/2.};
298     if(!(GetITSgeom()->IsShapeDefined(kSDD)))
299         GetITSgeom()->ReSetShape(kSDD, new AliITSgeomSDD256(3,(Float_t *)kDxyzSDD));
300     
301     const Float_t kDxyzSSD[] = {fgkSSDlengthSens/2, fgkSSDthickSens/2,fgkSSDwidthSens/2};
302     if(!(GetITSgeom()->IsShapeDefined(kSSD)))
303         GetITSgeom()->ReSetShape(kSSD,new AliITSgeomSSD75and275(3,(Float_t *)kDxyzSSD));
304     
305     // Creating the matrices in AliITSgeom for each sensitive volume
306     // (like in AliITSv11GeometrySDD) mln
307     // Here, each layer is one detector
308     
309     char layerName[30];
310     Int_t startMod = 0;
311     
312     // SPD
313     for (Int_t i=0; i<fNspd;i++) {
314         sprintf(layerName, "ITSspdWafer_%i",i+1);
315         TGeoNode *layNode = fITSmotherVolume->GetNode(layerName);
316         if (layNode) {
317             TGeoHMatrix layMatrix(*layNode->GetMatrix());           
318             Double_t *trans  = layMatrix.GetTranslation();
319             Double_t *r      = layMatrix.GetRotationMatrix();
320             Double_t rot[10] = {r[0],r[1],r[2],
321                                 r[3],r[4],r[5],
322                                 r[6],r[7],r[8], 1.0};
323             Int_t iDet = 1;
324             if ((i+1==2)||(i+1==4)) iDet = 2;
325             Int_t iLad = 1;
326             Int_t iLay = 1;
327             if (i+1>2) iLay = 2;
328             GetITSgeom()->CreatMatrix(startMod,iLay,iLad,iDet,kSPD,trans,rot);
329             startMod++;
330         };
331     };
332     
333     // SDD
334     for (Int_t i=0; i<fNsdd;i++) {
335         sprintf(layerName, "ITSsddWafer_%i",i+fNspd+1);
336         TGeoNode *layNode = fITSmotherVolume->GetNode(layerName);
337         if (layNode) {
338             TGeoHMatrix layMatrix(*layNode->GetMatrix());
339             Double_t *trans  = layMatrix.GetTranslation();
340             Double_t *r      = layMatrix.GetRotationMatrix();
341             Double_t rot[10] = {r[0],r[1],r[2],
342                                 r[3],r[4],r[5],
343                                 r[6],r[7],r[8], 1.0};
344             Int_t iDet = 1;
345             Int_t iLad = 1;
346             Int_t iLay = fNspd-1+i;
347             GetITSgeom()->CreatMatrix(startMod,iLay,iLad,iDet,kSDD,trans,rot);
348             startMod++;
349         };
350     };
351     
352     // SSD
353     for (Int_t i=0; i<fNssd;i++) {
354         sprintf(layerName, "ITSssdWafer_%i",i+fNspd+fNsdd+1);
355         TGeoNode *layNode = fITSmotherVolume->GetNode(layerName);
356         if (layNode) {
357             TGeoHMatrix layMatrix(*layNode->GetMatrix());           
358             Double_t *trans  = layMatrix.GetTranslation();
359             Double_t *r      = layMatrix.GetRotationMatrix();
360             Double_t rot[10] = {r[0],r[1],r[2],
361                                 r[3],r[4],r[5],
362                                 r[6],r[7],r[8], 1.0};
363             Int_t iDet = 1;
364             if ((i+1==2)||(i+1==4)) iDet = 2;
365             Int_t iLad = 1;
366             Int_t iLay = 5;
367             if (i+1>2) iLay = 6;
368             GetITSgeom()->CreatMatrix(startMod,iLay,iLad,iDet,kSSD,trans,rot);
369             startMod++;
370         };
371     };
372     
373     return;
374 }
375
376 //______________________________________________________________________
377 void AliITSvBeamTestITS04::SetDefaults()
378 {
379     // (from AliITSv11) mln
380     
381     const Float_t kconv = 1.0e+04; // convert cm to microns
382     
383     if(!fDetTypeSim) fDetTypeSim = new AliITSDetTypeSim();
384     fDetTypeSim->SetITSgeom(GetITSgeom());
385     fDetTypeSim->ResetCalibrationArray();
386     fDetTypeSim->ResetSegmentation();
387  
388     AliITSgeomSPD *s0;
389     AliITSgeomSDD *s1;
390     AliITSgeomSSD *s2;
391     Int_t i;
392     Float_t bx[256],bz[280];
393
394     // If fGeomDetIn is set true the geometry will
395     // be initialised from file (JC)
396     if(GetITSgeom()!=0) SetITSgeom(0x0);
397     AliITSgeom* geom = new AliITSgeom();
398     SetITSgeom(geom);
399     if(fGeomDetIn) GetITSgeom()->ReadNewFile(fRead);
400     if(!fGeomDetIn) this->InitAliITSgeom();
401     if(fGeomDetOut) GetITSgeom()->WriteNewFile(fWrite);
402
403    
404     // SPD
405
406     s0 = (AliITSgeomSPD*) GetITSgeom()->GetShape(kSPD);// Get shape info.
407     if (s0) {
408         AliITSCalibration *resp0=new AliITSCalibrationSPD();
409         SetCalibrationModel(kSPD,resp0);
410
411         AliITSsegmentationSPD *seg0=new AliITSsegmentationSPD(GetITSgeom());
412         seg0->SetDetSize(s0->GetDx()*2.*kconv, // base this on AliITSgeomSPD
413                          s0->GetDz()*2.*kconv, // for now.
414                          s0->GetDy()*2.*kconv);// x,z,y full width in microns.
415         seg0->SetNPads(256,160);               // Number of Bins in x and z
416         for(i=000;i<256;i++) bx[i] =  50.0;    // in x all are 50 microns.
417         for(i=000;i<160;i++) bz[i] = 425.0;    // most are 425 microns except below
418         for(i=160;i<280;i++) bz[i] =   0.0;    // Outside of detector.
419         bz[ 31] = bz[ 32] = 625.0;             // first chip boundry
420         bz[ 63] = bz[ 64] = 625.0;             // first chip boundry
421         bz[ 95] = bz[ 96] = 625.0;             // first chip boundry
422         bz[127] = bz[128] = 625.0;             // first chip boundry
423         bz[160] = 425.0;                       // Set so that there is no zero pixel size for fNz.
424         seg0->SetBinSize(bx,bz);               // Based on AliITSgeomSPD for now.
425         SetSegmentationModel(kSPD,seg0);
426         // set digit and raw cluster classes to be used
427         const char *kData0=(fDetTypeSim->GetCalibrationModel(kSPD))->DataType();
428         if (strstr(kData0,"real")) 
429           fDetTypeSim->SetDigitClassName(kSPD,"AliITSdigit");
430         else fDetTypeSim->SetDigitClassName(kSPD,"AliITSdigitSPD");
431     }
432   
433     // SDD
434    
435     s1 = (AliITSgeomSDD*) GetITSgeom()->GetShape(kSDD);// Get shape info.
436     if (s1) {
437       AliITSCalibrationSDD *resp1=new AliITSCalibrationSDD("simulated");
438       SetCalibrationModel(kSDD,resp1);
439       AliITSsegmentationSDD *seg1=new AliITSsegmentationSDD(GetITSgeom(),resp1);
440       seg1->SetDetSize(s1->GetDx()*kconv, // base this on AliITSgeomSDD
441                        s1->GetDz()*4.*kconv, // for now.
442                        s1->GetDy()*4.*kconv); // x,z,y full width in microns.
443       seg1->SetNPads(256,256);// Use AliITSgeomSDD for now
444       SetSegmentationModel(kSDD,seg1);
445       const char *kData1=(fDetTypeSim->GetCalibrationModel(kSDD))->DataType();
446       const char *kopt=resp1->GetZeroSuppOption();
447       if((!strstr(kopt,"2D")) && (!strstr(kopt,"1D")) || strstr(kData1,"real") ){
448         fDetTypeSim->SetDigitClassName(kSDD,"AliITSdigit");
449         } else fDetTypeSim->SetDigitClassName(kSDD,"AliITSdigitSDD");
450     }
451     
452     // SSD
453     
454     s2 = (AliITSgeomSSD*) GetITSgeom()->GetShape(kSSD);// Get shape info. Do it this way for now.
455     if (s2) {
456       AliITSCalibration *resp2=new AliITSCalibrationSSD("simulated");
457       SetCalibrationModel(kSSD,resp2);
458
459       AliITSsegmentationSSD *seg2=new AliITSsegmentationSSD(GetITSgeom());
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=(fDetTypeSim->GetCalibrationModel(kSSD))->DataType();
470       if(strstr(kData2,"real") ) fDetTypeSim->SetDigitClassName(kSSD,"AliITSdigit");
471       else fDetTypeSim->SetDigitClassName(kSSD,"AliITSdigitSSD");
472     }
473     
474   if(fgkNTYPES>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   }