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