]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSvSDD03.cxx
Modifications to reduce the size in memory (E.Crescio)
[u/mrichter/AliRoot.git] / ITS / AliITSvSDD03.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /*
17  $Id$ 
18 */
19 /////////////////////////////////////////////////////////////////
20 //  Class for the SDD beam test August2004                     //
21 //                                                             //
22 //                                                             //
23 /////////////////////////////////////////////////////////////////
24
25 #include <TGeometry.h>
26 #include <TNode.h>
27 #include <TBRIK.h>
28 #include <TLorentzVector.h>
29 #include <TVirtualMC.h>
30
31 #include "AliMC.h"
32 #include "AliRun.h"
33 #include "AliMagF.h"
34 #include "AliITSGeant3Geometry.h"
35 #include "AliTrackReference.h"
36 #include "AliITShit.h"
37 #include "AliITS.h"
38 #include "AliITSvSDD03.h"
39 #include "AliITSgeom.h"
40 #include "AliITSgeomSDD.h"
41 #include "AliITSgeomSSD.h"
42 #include "AliITSDetTypeSim.h"
43 #include "AliITSresponseSPD.h"
44 #include "AliITSresponseSDD.h"
45 #include "AliITSresponseSSD.h"
46 #include "AliITSsegmentationSPD.h"
47 #include "AliITSsegmentationSDD.h"
48 #include "AliITSsegmentationSSD.h"
49 #include "AliITSsimulationSPDdubna.h"
50 #include "AliITSsimulationSDD.h"
51 #include "AliITSsimulationSSD.h"
52
53
54 ClassImp(AliITSvSDD03)
55
56 //______________________________________________________________________
57 AliITSvSDD03::AliITSvSDD03() :
58 AliITS(),
59 fGeomDetOut(kFALSE),
60 fGeomDetIn(kFALSE),
61 fMajorVersion(1),
62 fMinorVersion(2),
63 fEuclidGeomDet(),
64 fRead(),
65 fWrite(),
66 fDet1(300.0),
67 fDet2(300.0),
68 fChip1(300.0),
69 fChip2(300.0),
70 fIDMother(0),
71 fYear(2003){
72     ////////////////////////////////////////////////////////////////////////
73     // Standard default constructor for the ITS SDD test beam 2002 version 1.
74     // Inputs:
75     //    none.
76     // Outputs:
77     //    none.
78     // Return:
79     //    A default created class.
80     ////////////////////////////////////////////////////////////////////////
81     Int_t i;
82
83     fIdN          = 0;
84     fIdName       = 0;
85     fIdSens       = 0;
86     fEuclidOut    = kFALSE; // Don't write Euclide file
87     for(i=0;i<60;i++) fRead[i] = '\0';
88     for(i=0;i<60;i++) fWrite[i] = '\0';
89     for(i=0;i<60;i++) fEuclidGeomDet[i] = '\0';
90     fTarg=kNoTarg;
91     fTargThick=0;
92 }
93 //______________________________________________________________________
94 AliITSvSDD03::AliITSvSDD03(const char *title,Int_t year):
95 AliITS("ITS", title),
96 fGeomDetOut(kFALSE),
97 fGeomDetIn(kFALSE),
98 fMajorVersion(1),
99 fMinorVersion(2),
100 fEuclidGeomDet(),
101 fRead(),
102 fWrite(),
103 fDet1(300.0),
104 fDet2(300.0),
105 fChip1(300.0),
106 fChip2(300.0),
107 fIDMother(0),
108 fYear(2003){
109     ////////////////////////////////////////////////////////////////////////
110     //    Standard constructor for the ITS SDD testbeam 2002 version 1.
111     // Inputs:
112     //    const char *title    title for this ITS geometry.
113     // Outputs:
114     //    none.
115     // Return:
116     //    A standard created class.
117     ////////////////////////////////////////////////////////////////////////
118     Int_t i;
119
120     fIdN = 3;
121     fIdName = new TString[fIdN];
122     fIdName[0] = "IMBS";
123     fIdName[1] = "ITST";
124     fIdName[2] = "ISNT";
125     fIdSens    = new Int_t[fIdN];
126     for(i=0;i<fIdN;i++) fIdSens[i] = 0;
127     fEuclidOut    = kFALSE; // Don't write Euclide file
128     fYear         = year;
129     fTarg=kNoTarg;
130     fTargThick=0;    
131     SetThicknessDet1();
132     SetThicknessDet2();
133     SetThicknessChip1();
134     SetThicknessChip2();                         
135
136     fEuclidGeometry="$ALICE_ROOT/ITS/ITSgeometry_vSDD032.euc";
137     strncpy(fEuclidGeomDet,"$ALICE_ROOT/ITS/ITSgeometry_vSDD032.det",60);
138     strncpy(fRead,fEuclidGeomDet,60);
139     strncpy(fWrite,fEuclidGeomDet,60);
140 }
141 //______________________________________________________________________
142 AliITSvSDD03::AliITSvSDD03(const AliITSvSDD03 &source) :  AliITS(source){
143     ////////////////////////////////////////////////////////////////////////
144     //     Copy Constructor for ITS SDD test beam 2002 version 1.
145     // This class is not to be copied. Function only dummy.
146     // Inputs:
147     //    const AliITSvSDD03 &source   The class to be copied
148     // Outputs:
149     //    none.
150     // Return:
151     //    A warning message.
152     ////////////////////////////////////////////////////////////////////////
153     if(&source == this) return;
154     Warning("Copy Constructor","Not allowed to copy AliITSvSDD03");
155     return;
156 }
157 //______________________________________________________________________
158 AliITSvSDD03& AliITSvSDD03::operator=(const AliITSvSDD03 &source){
159     ////////////////////////////////////////////////////////////////////////
160     //    Assignment operator for the ITS SDD test beam 2002 version 1.
161     // This class is not to be copied. Function only dummy.
162     // Inputs:
163     //    const AliITSvSDD03 &source   The class to be copied
164     // Outputs:
165     //    none.
166     // Return:
167     //    A Warning message
168     ////////////////////////////////////////////////////////////////////////
169     if(&source == this) return *this;
170     Warning("= operator","Not allowed to copy AliITSvSDD03");
171     return *this;
172 }
173 //______________________________________________________________________
174 AliITSvSDD03::~AliITSvSDD03() {
175     ////////////////////////////////////////////////////////////////////////
176     //    Standard destructor for the ITS SDD test beam 2002 version 1.
177     // Inputs:
178     //    none.
179     // Outputs:
180     //    none.
181     // Return:
182     //    none.
183     ////////////////////////////////////////////////////////////////////////
184 }
185 //______________________________________________________________________
186 void AliITSvSDD03::BuildGeometry(){
187     ////////////////////////////////////////////////////////////////////////
188     //    Geometry builder for the ITS SDD test beam 2002 version 1.
189     //    ALIC    ALICE Mother Volume
190     //     |- ITSV     ITS Mother Volume
191     //         |- IDET       Detector under Test
192     //         |   |- ITS0       SDD Si Chip
193     //         |   |  |- ITST      SDD Sensitivve Volume
194     //         |   |- IPC0 *5    Readout chip
195     //         |- ITEL *4    SDD Telescope
196     //             |- IMB0       SDD Si Chip
197     //             |   |- IMBS     SDD Sensitive volume
198     //             |- ICMB       Chip MiniBus.
199     // Inputs:
200     //    none.
201     // Outputs:
202     //    none.
203     // Return:
204     //    none.
205     ////////////////////////////////////////////////////////////////////////
206     // Get the top alice volume.
207     TNode *nALIC = gAlice->GetGeometry()->GetNode("alice");
208     nALIC->cd();
209
210     // Define ITS Mother Volume
211     Float_t data[3];
212     Float_t ddettest=200.0E-4,ddettelescope=300.0E-4;
213     //Float_t yposition= 0.0;
214     TRotMatrix *r0 = new TRotMatrix("ITSidrotm0","ITSidrotm0",
215                                     90.0,0,0.0,0,90.0,270.0);
216     data[0] = 10.0;
217     data[1] = 10.0;
218     data[2] = 100.0;
219     TBRIK *sITSVshape =new TBRIK("ITSVshape","ITS Logical Mother Volume","Air",
220                                  data[0],data[1],data[2]);
221     TNode *sITSV = new TNode("ITSV","ITS Mother Volume",sITSVshape,
222                             0.0,0.0,0.0,0,0);
223     sITSV->cd(); // set ourselve into ITSV subvolume of ALIC
224
225     // SSD part of telescope (MiniBuS)
226     data[0] = 1.06;
227     data[1] = 0.5*ddettelescope;
228     data[2] = 1.1;
229     TBRIK *sIMB0shape = new TBRIK("IMB0shape","SDD wafer","Si",
230                                  data[0],data[1],data[2]);
231     Float_t detMiniBusX,detMiniBusY,detMiniBusZ;
232     data[0] = detMiniBusX = 0.5*384*50.0E-4;
233     data[1] = detMiniBusY = 0.5*ddettelescope;
234     data[2] = detMiniBusZ = 1.0;
235     TBRIK *sIMBSshape = new TBRIK("IMBSshape","SDD Sensitive volume","Si",
236                                  data[0],data[1],data[2]);
237
238     data[0] = 1.36;
239     data[1] = 0.47;
240     data[2] = 1.36;
241     TBRIK *sITELshape = new TBRIK("ITELshape","ITELshape","Air",
242                                  data[0],data[1],data[2]);
243
244
245     // SDD under test
246     Float_t spdX,spdY,spdZ;
247     data[0] = 3.62500;
248     data[1] = 0.5*ddettest;
249     data[2] = 4.37940;
250     TBRIK *sITS0shape = new TBRIK("ITS0shape","SDD wafer","Si",
251                                  data[0],data[1],data[2]); // contains detector
252     data[0] = spdX = 3.50860;
253     data[1] = spdY = 0.5*ddettest;
254     data[2] = spdZ = 3.76320;
255     TBRIK *sITSTshape = new TBRIK("ITSTshape","SDD sensitive volume","Si",
256                                  data[0],data[1],data[2]);
257
258     data[0] = 4.2;
259     data[1] = 0.52;
260     data[2] = 5.2;
261     TBRIK *sIDETshape = new TBRIK("IDETshape","Detector Under Test","Air",
262                                  data[0],data[1],data[2]);
263
264
265     // Place volumes in geometry
266     char name[20],title[50];
267
268     //place SDD under test
269     Double_t px=0.0,py=0.0;
270     Double_t pz[2]={0.0,5.2};
271     TNode *nIDET[2],*nITS0[2],*nITST[2];
272     for(Int_t i=0;i<2;i++){
273         sITSV->cd();
274         sprintf(name,"IDET%d",i);
275         sprintf(title,"SDD #%d under test",i+1);
276         nIDET[i] = new TNode(name,title,sIDETshape,px,py,pz[i],r0,0);
277         nIDET[i]->cd();
278         nITS0[i] = new TNode("ITS0","SDD wafer",sITS0shape,0.0,0.0,0.0,0,0);
279         nITS0[i]->cd();
280         nITST[i] = new TNode("ITST","SDD sensitive volume",sITSTshape,
281                             0.0,0.0,0.0,0,0);
282         nITST[i]->SetLineColor(kYellow);
283         fNodes->Add(nITST[i]);
284     } // end for i
285
286     //place SSD telescope planes
287     Double_t qx=0.0,qy=0.0;
288     Double_t qz[10]={-58.4,-57.4,-50.4,-49.4,60.1,61.1,68.4,69.4,87.7,88.7};
289     TNode *nITEL[10],*nIMB0[10],*nIMBS[10];
290     for(Int_t i=0;i<10;i++){
291         sITSV->cd();
292         sprintf(name,"ITEL%d",i);
293         sprintf(title,"Test beam telescope element #%d",i+1);
294         nITEL[i] = new TNode(name,title,sITELshape,qx,qy,qz[i],r0,0);
295         nITEL[i]->cd();
296         nIMB0[i] = new TNode("IMB0","Chip MiniBus",sIMB0shape,
297                             0.0, 0.0,0.0,0,0);
298         nIMB0[i]->cd();
299         nIMBS[i] = new TNode("IMBS","IMBS",sIMBSshape,0.0,0.0,0.0,0,0);
300         nIMBS[i]->SetLineColor(kGreen);
301         fNodes->Add(nIMBS[i]);
302     } // end for i
303     nALIC->cd();
304     sITSV->Draw();
305 }
306 //______________________________________________________________________
307 Int_t AliITSvSDD03::DecodeDetector(Int_t id,Int_t cpy,Int_t &lay,
308                                    Int_t &det,Int_t &lad) const{
309     // Given the Geant id and copy volume number, returns the layer, ladder,
310     // and detector number, allong with the module number of the detector
311     // involved. Returns -1 and lay=0, lad=0, and det=0 if not a sensitive 
312     // volume.
313     // Inputs:
314     //    Int_t id    Geometry volume id number
315     //    Int_t cpy   Geometry copy number
316     // Outputs:
317     //    Int_t lay   ITS layer number
318     //    Int_t lad   ITS ladder number
319     //    Int_t det   ITS detector number
320     // Return:
321     //    Int_t module number.
322     Int_t mod;
323
324     lay = 0; lad = 0; det = 0; mod = -1;
325     if(id==fIdSens[0]){ // Volume name is IMBS (ITEL)
326         lad = 1; det = 1;
327         lay = cpy;
328         if(cpy>4) lay+=2;
329         mod = lay-1;
330         return mod;
331     }// end if
332     if(id==fIdSens[1]){ // Volume name is ITST (IDet)
333         lad = 1; det = 1;lay = cpy+4; 
334         mod = lay-1;
335         return mod;
336     }// end if
337     return mod;
338 }
339 //______________________________________________________________________
340 void AliITSvSDD03::CreateGeometry(){
341     ////////////////////////////////////////////////////////////////////////
342     //  This routine defines and Creates the geometry for version 1 of the ITS.
343     //    ALIC    ALICE Mother Volume
344     //     |- ITSV     ITS Mother Volume
345     //         |- IDET       Detector under Test (boxcontaining SDD)
346     //         |   |-IDAI        Air inside box
347     //         |       |- ITS0       SDD Si Chip
348     //         |          |- ITST      SDD Sensitivve Volume
349     //         |- ITEL *10   SSD Telescope (plastic box containting SSD's)
350     //         |   |- ITAI       Air inside box
351     //         |       |- IMB0       SDD Si Chip
352     //         |           |- IMBS     SDD Sensitive volume
353     //         |-ISNT*4    Sintilator triggers
354     // Inputs:
355     //    none.
356     // Outputs:
357     //    none.
358     // Return:
359     //    none.
360     ////////////////////////////////////////////////////////////////////////
361     Float_t data[49];
362     // Define media off-set
363     Int_t *idtmed = fIdtmed->GetArray()+1; // array of media indexes
364     Int_t idrotm[4]; // Array of rotation matrix indexes
365     //Float_t ddettest=200.0E-4,ddettelescope=300.0E-4;
366     //Float_t dchipMiniBus=750.0E-4,dchiptest=300.0E-4;
367     //Float_t yposition= 0.0;
368     const Float_t kmm=0.1,kcm=1.0,kmicm=kmm/1000.;
369     // Define Rotation-reflextion Matrixes needed
370     // 0 is the unit matrix
371     AliMatrix(idrotm[0], 90.0,0.0, 0.0,0.0, 90.0,270.0); // SDD and SSD X
372     AliMatrix(idrotm[1], 90.0,90.0, 0.0,180.0, 90.0,270.0); // SSD Y
373     AliMatrix(idrotm[2],90.0,90.0,90.0,180.0,0.0,0.0);  //Rotate about Z 90 degree
374
375     data[0] = 150.0*kmm;
376     data[1] = 150.0*kmm;
377     data[2] = 1100.0*kmm;
378     gMC->Gsvolu("ITSV","BOX ",idtmed[0],data,3);
379     gMC->Gspos("ITSV",1,"ALIC",0.0,0.0,0.0,0,"ONLY");
380
381
382     // Crossed sintilator triggers (2 in front 2 in back)
383     data[0] = 10.0*kcm;
384     data[1] = 2.0*kcm;
385     data[2] = 2.0*kmm;
386     gMC->Gsvolu("ISNT","BOX ",idtmed[2],data,3);
387     gMC->Gspos("ISNT",1,"ITSV",0.0,0.0,-950.0*kmm,0,"ONLY");
388     gMC->Gspos("ISNT",2,"ITSV",0.0,0.0,-950.0*kmm-data[2],idrotm[2],"ONLY");
389     gMC->Gspos("ISNT",3,"ITSV",0.0,0.0,950.0*kmm+data[2],0,"ONLY");
390     gMC->Gspos("ISNT",4,"ITSV",0.0,0.0,950.0*kmm,idrotm[2],"ONLY");
391
392
393 ////Create Volumes
394
395     // SSD part of telescope (MiniBuS)
396     Float_t detMiniBusX,detMiniBusY,detMiniBusZ;
397     data[0] = detMiniBusX = 10600.0*kmicm;
398     data[1] = detMiniBusY = 0.150*kmm;
399     data[2] = detMiniBusZ = 1.1*kcm;
400     gMC->Gsvolu("IMB0", "BOX ", idtmed[1], data, 3);   // contains detector
401     data[0] = 0.5*384*50*kmicm;
402     data[1] = 0.1499*kmm;
403     data[2] = 1.0*kcm;
404     gMC->Gsvolu("IMBS","BOX ",idtmed[1],data,3); // sensitive detector volume
405     gMC->Gspos("IMBS",1,"IMB0",0.0,0.0,0.0,0,"ONLY"); // place IMBS inside
406     // Box containing SSD's
407     data[0] = 11600.0*kmicm;
408     data[1] = 0.450*kcm;
409     data[2] = 1.16*kcm;
410     gMC->Gsvolu("ITAI","BOX ",idtmed[0],data,3);
411     // Plastic box size = insize + thickness.
412     data[0] = data[0] + 2.0*kmm;
413     data[1] = data[1] + 200.0*kmicm;
414     data[2] = data[2] + 2.0*kmm;
415     gMC->Gsvolu("ITEL","BOX ",idtmed[3],data,3);
416     gMC->Gspos("ITAI",1,"ITEL",0.0,0.0,0.0,0,"ONLY");
417     gMC->Gspos("IMB0",1,"ITAI",0.0,0.0,0.0,0,"ONLY");
418
419     // SDD under test
420     Float_t sddX,sddY,sddZ;
421     data[0] = sddX = 3.62500*kcm;
422     data[1] = sddY = 0.1500*kmm;
423     data[2] = sddZ = 4.37940*kcm;
424     gMC->Gsvolu("ITS0", "BOX ", idtmed[1], data, 3);   // contains detector
425     data[0] = 3.50860*kcm;
426     data[1] = 0.1499*kmm;
427     data[2] = 3.76320*kcm;
428     gMC->Gsvolu("ITST","BOX ",idtmed[1],data,3);// sensitive detecor volume
429     gMC->Gspos("ITST",1,"ITS0",0.0,0.0,0.0,0,"ONLY"); // place ITST inside
430     // Box containing SDD under test
431     data[0] = 4.0*kcm;
432     data[1] = 0.5*kcm;
433     data[2] = 5.0*kcm;
434     gMC->Gsvolu("IDAI","BOX ",idtmed[0],data,3);
435     data[0] = data[0] + 2.0*kmm;
436     data[1] = data[1] + 200.0*kmicm;
437     data[2] = data[2] + 2.0*kmm;
438     gMC->Gsvolu("IDET","BOX ",idtmed[3],data,3);
439     gMC->Gspos("IDAI",1,"IDET",0.0,0.0,0.0,0,"ONLY");
440     gMC->Gspos("ITS0",1,"IDAI",0.0,0.0,0.0,0,"ONLY");
441
442
443 //// Position detectors, Beam Axis Z, X to the right, Y up to the sky.
444     // Upsteram planes of the telescope    
445     Float_t p00X,p00Y,p00Z,p01X,p01Y,p01Z,p10X,p10Y,p10Z,p11X,p11Y,p11Z;
446     p00X = 0.0*kcm;
447     p00Y = 0.0*kcm;
448     p00Z = -584*kmm;
449     gMC->Gspos("ITEL",1,"ITSV",p00X,p00Y,p00Z,idrotm[0],"ONLY");//SSD X
450     p01X = 0.0*kcm;
451     p01Y = 0.0*kcm;
452     p01Z = -574*kmm;
453     gMC->Gspos("ITEL",2,"ITSV",p01X,p01Y,p01Z,idrotm[1],"ONLY");//SSD Y
454     p01X = 0.0*kcm;
455     p01Y = 0.0*kcm;
456     p01Z = -504*kmm;
457     gMC->Gspos("ITEL",3,"ITSV",p01X,p01Y,p01Z,idrotm[0],"ONLY");//SSD X
458     p01X = 0.0*kcm;
459     p01Y = 0.0*kcm;
460     p01Z = -494*kmm;
461     gMC->Gspos("ITEL",4,"ITSV",p01X,p01Y,p01Z,idrotm[1],"ONLY");//SSD Y
462
463     // Downstream planes of the telescope
464     p10X = 0.0*kcm;
465     p10Y = 0.0*kcm;
466     p10Z = +601.0*kmm; 
467     gMC->Gspos("ITEL",5,"ITSV",p10X,p10Y,p10Z,idrotm[0],"ONLY");//SSD X
468     p11X = 0.0*kcm;
469     p11Y = 0.0*kcm;
470     p11Z = +610.0*kmm; //611.0
471     gMC->Gspos("ITEL",6,"ITSV",p11X,p11Y,p11Z,idrotm[1],"ONLY");//SSD Y
472     p11X = 0.0*kcm;
473     p11Y = 0.0*kcm;
474     p11Z = +684.0*kmm;
475     gMC->Gspos("ITEL",7,"ITSV",p11X,p11Y,p11Z,idrotm[0],"ONLY");//SSD X
476     p11X = 0.0*kcm;
477     p11Y = 0.0*kcm;
478     p11Z = +694.0*kmm;
479     gMC->Gspos("ITEL",8,"ITSV",p11X,p11Y,p11Z,idrotm[1],"ONLY");//SSD Y
480     p11X = 0.0*kcm;
481     p11Y = 0.0*kcm;
482     p11Z = +877.0*kmm;
483     gMC->Gspos("ITEL",9,"ITSV",p11X,p11Y,p11Z,idrotm[0],"ONLY");//SSD X
484     p11X = 0.0*kcm;
485     p11Y = 0.0*kcm;
486     p11Z = +887.0*kmm;
487     gMC->Gspos("ITEL",10,"ITSV",p11X,p11Y,p11Z,idrotm[1],"ONLY");//SSD Y
488
489     // SDDs 
490     Float_t pdet1X,pdet1Y,pdet1Z;
491     Float_t pdet2X,pdet2Y,pdet2Z;
492     pdet1X = 0.0*kcm;
493     pdet1Y = 0.0*kcm;
494     pdet1Z = 0.0*kcm;
495     gMC->Gspos("IDET",1,"ITSV",pdet1X,pdet1Y,pdet1Z,idrotm[0],"ONLY");// Detector1
496     pdet2X = 0.0*kcm;
497     pdet2Y = 0.0*kcm;
498     pdet2Z = 52*kmm; //52
499     gMC->Gspos("IDET",2,"ITSV",pdet2X,pdet2Y,pdet2Z,idrotm[0],"ONLY");// Detector2
500
501 // Target definition and placement
502     if(fTarg){
503       data[0] = 30*kmm;
504       data[1] = fTargThick*kmm;  // Target thickness
505       data[2] = 30*kmm;
506       gMC->Gsvolu("ITGT","BOX ",idtmed[fTarg],data,3);
507
508       Float_t a,z,dens,radl,absl;
509       Float_t* ubuf=0; Int_t nbuf;
510       char* ssss="";
511       gMC->Gfmate(idtmed[fTarg],ssss,a,z,dens,radl,absl,ubuf,nbuf);
512
513       Info("CreateGeometry","Target A=%f,  Z=%f,  dens=%f",a,z,dens);
514       Info("Creategeometry","Target thickness=%f mm",fTargThick);
515
516       Float_t ptgtX,ptgtY,ptgtZ;
517       ptgtX = 0.0*kcm;
518       ptgtY = 0.0*kcm;
519       ptgtZ = -50*kmm;
520       gMC->Gspos("ITGT",1,"ITSV",ptgtX,ptgtY,ptgtZ,idrotm[0],"ONLY");// Target
521     }else{
522       Info("CreateGeometry","No target defined");
523     }
524 }
525 //______________________________________________________________________
526 void AliITSvSDD03::CreateMaterials(){
527     ////////////////////////////////////////////////////////////////////////
528     //
529     // Create ITS SDD test beam materials
530     //     This function defines the default materials used in the Geant
531     // Monte Carlo simulations for the geometries AliITSv1, AliITSv3,
532     // AliITSvSDD03.
533     // In general it is automatically replaced by
534     // the CreatMaterials routine defined in AliITSv?. Should the function
535     // CreateMaterials not exist for the geometry version you are using this
536     // one is used. See the definition found in AliITSv5 or the other routine
537     // for a complete definition.
538     //
539     // Inputs:
540     //    none.
541     // Outputs:
542     //    none.
543     // Return:
544     //    none.
545     /////////////////////////////////////////////////////////////////////////
546     Float_t tmaxfdSi = 0.1; // Degree
547     Float_t stemaxSi = 0.0075; // cm //0.0075
548     Float_t deemaxSi = 0.1; // Fraction of particle's energy 0<deemax<=1
549     Float_t epsilSi  = 1.0E-4;//
550     Float_t stminSi  = 0.0; // cm "Default value used"
551
552     Float_t tmaxfdAir = 0.1; // Degree
553     Float_t stemaxAir = .10000E+01; // 1 cm  //cm
554     Float_t deemaxAir = 0.1; // Fraction of particle's energy 0<deemax<=1
555     Float_t epsilAir  = 1.0E-4;//
556     Float_t stminAir  = 0.0; // cm "Default value used"
557     Int_t   ifield = gAlice->Field()->Integ();
558     Float_t fieldm = gAlice->Field()->Max();
559     //
560
561     // AIR
562     Float_t aAir[4]={12.0107,14.0067,15.9994,39.948};
563     Float_t zAir[4]={6.,7.,8.,18.};
564     Float_t wAir[4]={0.000124,0.755267,0.231781,0.012827};
565     Float_t dAir = 1.20479E-3;
566     // Lucite/Plexiglass
567     Float_t aLuc[3] = {1.,12.,16.};
568     Float_t zLuc[3] = {1.,6.,8.};
569     Float_t wLuc[3] = {8.,5.,2.};
570     Float_t dLuc = 1.19;
571     // stainless steel
572     Float_t asteel[4] = { 55.847,51.9961,58.6934,28.0855 };
573     Float_t zsteel[4] = { 26.,24.,28.,14. };
574     Float_t wsteel[4] = { .715,.18,.1,.005 };
575     Float_t dsteel = 7.88;
576
577     AliMixture(1, "AIR$",aAir,zAir,dAir,4,wAir);
578     AliMaterial(2,"SI$",28.086,14.0,2.3300,9.3600,999.00);
579     AliMixture(3,"Sintilator$",aLuc,zLuc,dLuc,-3,wLuc);
580     AliMixture(4,"PlasticBox$",aLuc,zLuc,dLuc,-3,wLuc);
581     AliMaterial(5, "IRON$", 55.85, 26., 7.87, 1.76, 999.00);
582     AliMaterial(6, "LEAD$", 207.19, 82., 11.35, .56, 999.00);
583     AliMixture(7, "STAINLESS STEEL$", asteel, zsteel,dsteel, 4, wsteel);
584     AliMaterial(9, "C$", 12.011, 6., 2.265, 18.8, 999.00);
585     AliMaterial(10, "Al$", 26.98, 13., 2.70, 8.9, 999.00);
586     AliMaterial(11, "Be$", 9.012, 4., 1.848, 35.3, 999.00);
587     AliMaterial(12, "Ti$", 47.88, 22., 4.54, 3.56, 999.00);
588     AliMaterial(13, "Sn$", 118.69, 50., 7.31, 1.21, 999.00); 
589     AliMaterial(14, "Cu$", 63.55, 29., 8.96, 1.43, 999.00);
590     AliMaterial(15, "Ge$", 72.59, 32., 5.323, 2.30, 999.00);
591     AliMaterial(20, "W$", 183.85, 74., 19.3, 0.35, 999.00);
592  
593     AliMedium(1,"AIR$",1,0,ifield,fieldm,tmaxfdAir,stemaxAir,deemaxAir,
594               epsilAir,stminAir);
595     AliMedium(2,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
596               epsilSi,stminSi);
597     AliMedium(3,"Scintillator$",3,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
598               epsilSi,stminSi);
599     AliMedium(4,"PlasticBox$",4,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
600               epsilSi,stminSi);
601     AliMedium(5,"IRON$",5,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
602               epsilSi,stminSi);
603     AliMedium(6,"LEAD$",6,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
604               epsilSi,stminSi);
605     AliMedium(7,"StainlessSteel$",7,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
606               epsilSi,stminSi);
607
608     AliMedium(9,"C$",9,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,epsilSi,stminSi);
609     AliMedium(10,"Al$",10,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,epsilSi,stminSi);
610     AliMedium(11,"Be$",11,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,epsilSi,stminSi);
611     AliMedium(12,"Ti$",12,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,epsilSi,stminSi);
612     AliMedium(13,"Sn$",13,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,epsilSi,stminSi);
613     AliMedium(14,"Cu$",14,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,epsilSi,stminSi);
614     AliMedium(15,"Ge$",15,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,epsilSi,stminSi);
615     AliMedium(20,"W$",20,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,epsilSi,stminSi);
616     //dummy materials to avoid warning during simulation (galice.cuts)
617
618    AliMedium(21,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
619               epsilSi,stminSi);
620    AliMedium(25,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
621               epsilSi,stminSi);
622    AliMedium(26,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
623               epsilSi,stminSi);
624    AliMedium(27,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
625               epsilSi,stminSi);
626    AliMedium(51,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
627               epsilSi,stminSi);
628    AliMedium(52,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
629               epsilSi,stminSi);
630    AliMedium(53,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
631               epsilSi,stminSi);
632    AliMedium(54,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
633               epsilSi,stminSi);
634    AliMedium(55,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
635               epsilSi,stminSi);
636    AliMedium(56,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
637               epsilSi,stminSi);
638    AliMedium(61,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
639               epsilSi,stminSi);
640    AliMedium(62,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
641               epsilSi,stminSi);
642    AliMedium(63,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
643               epsilSi,stminSi);
644    AliMedium(64,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
645               epsilSi,stminSi);
646    AliMedium(65,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
647               epsilSi,stminSi);
648    AliMedium(68,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
649               epsilSi,stminSi);
650    AliMedium(69,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
651               epsilSi,stminSi);
652    AliMedium(70,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
653               epsilSi,stminSi);
654    AliMedium(71,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
655               epsilSi,stminSi);
656    AliMedium(72,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
657               epsilSi,stminSi);
658    AliMedium(73,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
659               epsilSi,stminSi);
660    AliMedium(74,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
661               epsilSi,stminSi);
662    AliMedium(75,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
663               epsilSi,stminSi);
664    AliMedium(76,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
665               epsilSi,stminSi);
666    AliMedium(77,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
667               epsilSi,stminSi);
668    AliMedium(78,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
669               epsilSi,stminSi);
670    AliMedium(79,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
671               epsilSi,stminSi);
672    AliMedium(80,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
673               epsilSi,stminSi);
674    AliMedium(81,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
675               epsilSi,stminSi);
676    AliMedium(82,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
677               epsilSi,stminSi);
678    AliMedium(83,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
679               epsilSi,stminSi);
680    AliMedium(84,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
681               epsilSi,stminSi);
682    AliMedium(85,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
683               epsilSi,stminSi);
684    AliMedium(90,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
685               epsilSi,stminSi);
686    AliMedium(91,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
687               epsilSi,stminSi);
688    AliMedium(92,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
689               epsilSi,stminSi);
690    AliMedium(93,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
691               epsilSi,stminSi);
692    AliMedium(94,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
693               epsilSi,stminSi);
694    AliMedium(95,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
695               epsilSi,stminSi);
696    AliMedium(96,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
697               epsilSi,stminSi);
698
699 }            
700 //______________________________________________________________________
701 void AliITSvSDD03::InitAliITSgeom(){
702     //     Based on the geometry tree defined in Geant 3.21, this
703     // routine initilizes the Class AliITSgeom from the Geant 3.21 ITS geometry
704     // sturture.
705     // Inputs:
706     //    none.
707     // Outputs:
708     //    none.
709     // Return:
710     //    none.
711
712   
713   if(strcmp(gMC->GetName(),"TGeant3")) {
714     Error("InitAliITSgeom",
715           "Wrong Monte Carlo. InitAliITSgeom uses TGeant3 calls");
716     return;
717   } // end if
718   Info("InitAliITSgeom","Reading geometry transformation directly from Geant 3");
719   const Int_t knp=384;
720   const Float_t kpitch=50.E-4;/*cm*/
721   Float_t box[3]={0.5*kpitch*(Float_t)knp,150.E-4,1.0},p[knp+1],n[knp+1];
722   const Int_t kltypess = 2;
723   const Int_t knlayers = 12;
724   const Int_t kndeep = 6;
725   Int_t itsGeomTreeNames[kltypess][kndeep],lnam[20],lnum[20];
726   Int_t nlad[knlayers],ndet[knlayers];
727   Double_t t[3],r[10];
728   Float_t  par[20],att[20];
729   Int_t    npar,natt,idshape,imat,imed,id;
730   AliITSGeant3Geometry *ig = new AliITSGeant3Geometry();
731   Int_t mod,typ,lay,lad,det,cpy,i,j,k;
732   Char_t names[kltypess][kndeep][4];
733   Int_t itsGeomTreeCopys[kltypess][kndeep];
734   Char_t *namesA[kltypess][kndeep] = {
735     {"ALIC","ITSV","ITEL","ITAI","IMB0","IMBS"}, 
736     {"ALIC","ITSV","IDET","IDAI","ITS0","ITST"}};
737   Int_t itsGeomTreeCopysA[kltypess][kndeep]= {{1,1,10,1,1,1},
738                                             {1,1,2,1,1,1}};
739   for(i=0;i<kltypess;i++)for(j=0;j<kndeep;j++){
740     for(k=0;k<4;k++) names[i][j][k] = namesA[i][j][k];
741     itsGeomTreeCopys[i][j] = itsGeomTreeCopysA[i][j];
742   } // end for i,j
743   p[0]=-box[0];
744   n[0]=box[0];
745   for(i=1;i<knp;i++){// Fill in anode and cathode strip locations (lower edge)
746     p[i] =p[i-1]+kpitch;
747     n[i] =n[i-1]-kpitch;
748   } // end for i
749   p[knp+1]=box[0];
750   n[knp+1]=-box[0];
751
752   for(i=0;i<20;i++) lnam[i] = lnum[i] = 0;
753   for(i=0;i<kltypess;i++)for(j=0;j<kndeep;j++) 
754     strncpy((char*) &itsGeomTreeNames[i][j],names[i][j],4);
755   //    itsGeomTreeNames[i][j] = ig->StringToInt(names[i][j]);
756   mod = knlayers;
757   if(GetITSgeom()!=0) SetITSgeom(0x0);
758   for(Int_t iMod=0; iMod<knlayers; iMod++){
759     nlad[iMod]=1; ndet[iMod]=1;
760   }
761   
762   AliITSgeom* geom = new AliITSgeom(0,knlayers,nlad,ndet,mod);
763   SetITSgeom(geom);
764   fIdSens[0] = 0; fIdSens[1] = 1; // Properly reset in Init later.
765   for(typ=1;typ<=kltypess;typ++){
766     for(j=0;j<kndeep;j++) lnam[j] = itsGeomTreeNames[typ-1][j];
767     for(j=0;j<kndeep;j++) lnum[j] = itsGeomTreeCopys[typ-1][j];
768     if(typ == 1) id = fIdSens[0];
769     else id = fIdSens[1];
770     lad = 1;
771     det = 1;
772     for(cpy=1;cpy<=itsGeomTreeCopys[typ-1][2];cpy++){
773       mod = DecodeDetector(id,cpy,lay,lad,det);
774       //lunm[2] is the copy number, from 1 to 10 for SSD
775       // 1,2 for SDD, depending on the module number
776       if(mod<=3) lnum[2]=mod+1;
777       if(mod>=6) lnum[2]=mod-1;
778       if(mod==4) lnum[2]=1;
779       if(mod==5) lnum[2]=2;
780       ig->GetGeometry(kndeep,lnam,lnum,t,r,idshape,npar,natt,par,att,
781                       imat,imed);
782
783       switch (typ){
784       case 2:
785         GetITSgeom()->CreatMatrix(mod,lay,lad,det,kSDD,t,r);
786         if(!(GetITSgeom()->IsShapeDefined((Int_t)kSDD))){
787           GetITSgeom()->ReSetShape(kSDD,new AliITSgeomSDD256(npar,par));
788         } // end if
789         break;
790       case 1:
791         GetITSgeom()->CreatMatrix(mod,lay,lad,det,kSSD,t,r);
792         if(!(GetITSgeom()->IsShapeDefined((Int_t)kSSD))){
793           GetITSgeom()->ReSetShape(kSSD,new AliITSgeomSSD(box,0.0,0.0,
794                                                       knp+1,p,knp+1,n));
795         } // end if
796         break;
797       } // end switch
798     } // end for cpy
799   } // end for typ
800   return;
801 }
802 //______________________________________________________________________
803 void AliITSvSDD03::Init(){
804     ////////////////////////////////////////////////////////////////////////
805     //     Initialise the ITS after it has been created.
806     // Inputs:
807     //    none.
808     // Outputs:
809     //    none.
810     // Return:
811     //    none.
812     ////////////////////////////////////////////////////////////////////////
813
814
815     Info("Init","**********AliITSvSDD03 %d _Init *************",fMinorVersion);
816
817     if(fRead[0]=='\0') strncpy(fRead,fEuclidGeomDet,60);
818     if(fWrite[0]=='\0') strncpy(fWrite,fEuclidGeomDet,60);
819     if(GetITSgeom()!=0) SetITSgeom(0x0);
820     AliITSgeom* geom = new AliITSgeom();
821     SetITSgeom(geom);
822     if(fGeomDetIn) GetITSgeom()->ReadNewFile(fRead);
823     if(!fGeomDetIn) this->InitAliITSgeom();
824     if(fGeomDetOut) GetITSgeom()->WriteNewFile(fWrite);
825     AliITS::Init();
826     fIDMother = gMC->VolId("ITSV"); // ITS Mother Volume ID.
827
828 }
829 //______________________________________________________________________
830 void AliITSvSDD03::SetDefaults(){
831     // sets the default segmentation, response, digit and raw cluster classes
832     // Inputs:
833     //    none.
834     // Outputs:
835     //    none.
836     // Return:
837     //    none.
838
839     const Float_t kconv = 1.0e+04; // convert cm to microns
840
841     if(!fDetTypeSim) fDetTypeSim = new AliITSDetTypeSim();
842     fDetTypeSim->SetITSgeom(GetITSgeom());
843     fDetTypeSim->ResetResponse();
844     fDetTypeSim->ResetSegmentation();
845  
846     AliITSgeomSDD *s1;
847     AliITSgeomSSD *s2;
848     SetResponseModel(kSPD,new AliITSresponseSPD());
849     SetSegmentationModel(kSPD,new AliITSsegmentationSPD());
850     fDetTypeSim->SetDigitClassName(kSPD,"AliITSdigitSPD");
851
852     // SDD
853     s1 = (AliITSgeomSDD*) GetITSgeom()->GetShape(kSDD);// Get shape info. Do it this way for now.
854     AliITSresponseSDD *resp1=new AliITSresponseSDD("simulated");
855     SetResponseModel(kSDD,resp1);
856
857     AliITSsegmentationSDD *seg1=new AliITSsegmentationSDD(GetITSgeom(),resp1);
858     seg1->SetDetSize(s1->GetDx()*kconv, // base this on AliITSgeomSDD
859                      s1->GetDz()*2.*kconv, // for now.
860                      s1->GetDy()*2.*kconv); // x,z,y full width in microns.
861
862     seg1->SetNPads(256,256);// Use AliITSgeomSDD for now
863     SetSegmentationModel(kSDD,seg1);
864     const char *kData1=(fDetTypeSim->GetResponseModel(kSDD))->DataType();
865     const char *kopt=fDetTypeSim->GetResponseModel(kSDD)->ZeroSuppOption();
866     if((!strstr(kopt,"2D")) && (!strstr(kopt,"1D")) || strstr(kData1,"real") ){
867         fDetTypeSim->SetDigitClassName(kSDD,"AliITSdigit");
868     } else fDetTypeSim->SetDigitClassName(kSDD,"AliITSdigitSDD");
869
870     // SSD  Layer 5
871
872     s2 = (AliITSgeomSSD*) GetITSgeom()->GetShape(kSSD);// Get shape info. Do it this way for now.
873    
874     AliITSresponse *resp2= new AliITSresponseSSD("simulated");
875     SetResponseModel(kSSD,resp2);
876
877     AliITSsegmentationSSD *seg2=new AliITSsegmentationSSD(GetITSgeom());
878     seg2->SetDetSize(s2->GetDx()*2.*kconv, // base this on AliITSgeomSSD
879                      s2->GetDz()*2.*kconv, // for now.
880                      s2->GetDy()*2.*kconv); // x,z,y full width in microns.
881     seg2->SetPadSize(50.,0.); // strip x pitch in microns
882     seg2->SetNPads(384,0); // number of strips on each side.
883     seg2->SetLayer(5);
884     seg2->SetAngles(0.,0.); // strip angles rad P and N side.
885     seg2->SetAnglesLay5(0.,0.); // strip angles rad P and N side.
886     seg2->SetAnglesLay6(0.,0.); // strip angles rad P and N side.
887
888     SetSegmentationModel(kSSD,seg2); 
889     const char *kData2=(fDetTypeSim->GetResponseModel(kSSD))->DataType();
890     if(strstr(kData2,"real") ) fDetTypeSim->SetDigitClassName(kSSD,"AliITSdigit");
891     else fDetTypeSim->SetDigitClassName(kSSD,"AliITSdigitSSD");
892
893     if(fgkNTYPES>3){
894         Warning("SetDefaults",
895                 "Only the four basic detector types are initialised!");
896     }// end if
897     return;
898 }
899 //______________________________________________________________________
900 void AliITSvSDD03::SetDefaultSimulation(){
901     // sets the default simulation.
902     // Inputs:
903     //      none.
904     // Outputs:
905     //      none.
906     // Return:
907     //      none.
908
909   if(!fDetTypeSim) fDetTypeSim = new AliITSDetTypeSim();
910   AliITSsimulation *sim;
911   AliITSsegmentation *seg;
912   AliITSresponse *res;
913   //SPD
914   if(fDetTypeSim){
915     sim = fDetTypeSim->GetSimulationModel(kSPD);
916     if (!sim) {
917       seg =(AliITSsegmentation*)fDetTypeSim->GetSegmentationModel(kSPD);
918       if(seg==0) seg = new AliITSsegmentationSPD();
919       res = (AliITSresponse*)fDetTypeSim->GetResponseModel(kSPD);
920       if(res==0) res = new AliITSresponseSPD();
921       sim = new AliITSsimulationSPDdubna(seg,res,0);
922       SetSimulationModel(kSPD,sim);
923     }else{ // simulation exists, make sure it is set up properly.
924       sim->SetSegmentationModel((AliITSsegmentation*)fDetTypeSim->GetSegmentationModel(kSPD));
925       sim->SetResponseModel((AliITSresponse*)fDetTypeSim->GetResponseModel(kSPD));
926       sim->Init();
927     } // end if
928   } // end if iDetType
929   //SDD
930   if(fDetTypeSim){
931     sim = fDetTypeSim->GetSimulationModel(kSDD);
932     if (!sim) {
933       seg = (AliITSsegmentation*)fDetTypeSim->GetSegmentationModel(kSDD);
934       res = (AliITSresponse*)fDetTypeSim->GetResponseModel(kSDD);
935       sim = new AliITSsimulationSDD(seg,res);
936       SetSimulationModel(kSDD,sim);
937     }else{ // simulation exists, make sure it is set up properly.
938       sim->SetSegmentationModel((AliITSsegmentation*)fDetTypeSim->GetSegmentationModel(kSDD));
939       sim->SetResponseModel((AliITSresponse*)fDetTypeSim->GetResponseModel(kSDD));
940       
941       sim->Init();
942     } //end if
943   } // end if iDetType
944   //SSD
945   if(fDetTypeSim){
946     sim = fDetTypeSim->GetSimulationModel(kSSD);
947     if (!sim) {
948       seg = (AliITSsegmentation*)fDetTypeSim->GetSegmentationModel(kSSD);
949       res = (AliITSresponse*)fDetTypeSim->GetResponseModel(kSSD);
950       sim = new AliITSsimulationSSD(seg,res);
951       SetSimulationModel(kSSD,sim);
952     }else{ // simulation exists, make sure it is set up properly.
953       sim->SetSegmentationModel((AliITSsegmentation*)fDetTypeSim->GetSegmentationModel(kSSD));
954       sim->SetResponseModel((AliITSresponse*)fDetTypeSim->GetResponseModel(kSSD));
955       sim->Init();
956     } // end if
957   } // end if iDetType
958 }
959 //______________________________________________________________________
960 void AliITSvSDD03::DrawModule() const{
961     ////////////////////////////////////////////////////////////////////////
962     //     Draw a shaded view of the ITS SDD test beam version 1.
963     // Inputs:
964     //    none.
965     // Outputs:
966     //    none.
967     // Return:
968     //    none.
969     ////////////////////////////////////////////////////////////////////////
970     // Set everything unseen
971     gMC->Gsatt("*", "seen", -1);
972     // Set ALIC mother visible
973     gMC->Gsatt("ALIC","SEEN",0);
974     // Set ALIC ITS visible
975     gMC->Gsatt("ITSV","SEEN",1);
976     // Set ALIC Telescopes visible
977     gMC->Gsatt("ITEL","SEEN",1);
978     gMC->Gsatt("ITEL","colo",2);
979     // Set ALIC detetcor visible
980     gMC->Gsatt("IDET","SEEN",1);
981     gMC->Gsatt("IDET","colo",4);
982     // Set ALIC Scintillator visible
983     gMC->Gsatt("ISNT","SEEN",1);
984     gMC->Gsatt("ISNT","colo",3);
985     // Set Detector mother visible and drawn
986 //    gMC->Gsatt("ITS0","SEEN",1);
987     // Set minibus mother visible and drawn
988 //    gMC->Gsatt("IMB0","SEEN",1);
989
990     // Draw
991     gMC->Gdraw("alic", 60, 30, 180, 10,10, .12, .12);
992 }
993 //______________________________________________________________________
994 void AliITSvSDD03::StepManager(){
995     ////////////////////////////////////////////////////////////////////////
996     //    Called for every step in the ITS SDD test beam, then calles the 
997     // AliITShit class  creator with the information to be recoreded about
998     //  that hit.
999     //     The value of the macro ALIITSPRINTGEOM if set to 1 will allow the
1000     // printing of information to a file which can be used to create a .det
1001     // file read in by the routine CreateGeometry(). If set to 0 or any other
1002     // value except 1, the default behavior, then no such file is created nor
1003     // it the extra variables and the like used in the printing allocated.
1004     // Inputs:
1005     //    none.
1006     // Outputs:
1007     //    none.
1008     // Return:
1009     //    none.
1010     ////////////////////////////////////////////////////////////////////////
1011     Int_t  copy, id;
1012     static TLorentzVector position0;
1013     static Int_t stat0=0;
1014     if((id=gMC->CurrentVolID(copy) == fIDMother)&&
1015        (gMC->IsTrackEntering()||gMC->IsTrackExiting())){
1016         copy = fTrackReferences->GetEntriesFast();
1017         TClonesArray &lTR = *fTrackReferences;
1018         // Fill TrackReference structure with this new TrackReference.
1019         new(lTR[copy]) AliTrackReference(gAlice->GetMCApp()->
1020                                          GetCurrentTrackNumber());
1021     } // if Outer ITS mother Volume
1022     //if(!(this->IsActive())) return;
1023     if(!(gMC->TrackCharge())) return;
1024     Int_t   vol[5],copy3;
1025     TLorentzVector position, momentum;
1026     TClonesArray &lhits = *fHits;
1027     //
1028     // Fill hit structure.
1029     gMC->TrackPosition(position);
1030     gMC->TrackMomentum(momentum);
1031     id   = gMC->CurrentVolID(copy);
1032     if(id==fIdSens[0] || id==fIdSens[1]){ // Volumes "ITST" or "IMBS"
1033         copy = gMC->CurrentVolOffID(3,copy3);
1034         copy = DecodeDetector(id,copy3,vol[0],vol[1],vol[2]);
1035         //    cout << "0: mod,lay,lad,det="<<copy<<","<<vol[0]<<","<<vol[1]
1036         //     <<","<<vol[2]<<" name="<<gMC->CurrentVolName()
1037         //     <<","<<fIdSens[0]<<","<<fIdSens[1]<<","<<id<<" z="
1038         //     <<position.Z()<<endl;
1039     }else if(id==fIdSens[2]){ // "ISNT" Sintilator
1040         //cout << "1: id,copy="<<id<<","<<copy
1041         //     <<" name="<<gMC->CurrentVolName()<<" z="
1042         //     <<position.Z()<<endl;
1043         return; // Do nothing for now.
1044     }else{
1045         //cout << "2: id,copy="<<id<<","<<copy
1046         //     <<" name="<<gMC->CurrentVolName()<<" z="
1047         //     <<position.Z()<<endl;
1048         return;
1049     } // end if
1050     //
1051     // Track status
1052     vol[3] = 0;
1053     vol[4] = 0;
1054     if(gMC->IsTrackInside())      vol[3] +=  1;
1055     if(gMC->IsTrackEntering())    vol[3] +=  2;
1056     if(gMC->IsTrackExiting())     vol[3] +=  4;
1057     if(gMC->IsTrackOut())         vol[3] +=  8;
1058     if(gMC->IsTrackDisappeared()) vol[3] += 16;
1059     if(gMC->IsTrackStop())        vol[3] += 32;
1060     if(gMC->IsTrackAlive())       vol[3] += 64;
1061     //
1062     vol[4] = stat0;
1063     if(gMC->IsTrackEntering()){
1064         position0 = position;
1065         stat0 = vol[3];
1066         return;
1067     }else{
1068         new(lhits[fNhits++]) AliITShit(fIshunt,
1069                                    gAlice->GetMCApp()->GetCurrentTrackNumber(),
1070                                        vol,gMC->Edep(),gMC->TrackTime(),
1071                                        position,position0,momentum);
1072     }// end if gMC->IsTrackEntering()
1073     position0 = position;
1074     stat0 = vol[3];
1075
1076     return;
1077 }
1078