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