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