]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSvSDD03.cxx
added missing #include "AliMC.h" .
[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 #include <Riostream.h>
21 #include <stdio.h>
22 #include <stdlib.h>
23 #include <TMath.h>
24 #include <TGeometry.h>
25 #include <TNode.h>
26 #include <TTUBE.h>
27 #include <TTUBS.h>
28 #include <TPCON.h>
29 #include <TFile.h>    // only required for Tracking function?
30 #include <TCanvas.h>
31 #include <TObjArray.h>
32 #include <TLorentzVector.h>
33 #include <TObjString.h>
34 #include <TClonesArray.h>
35 #include <TBRIK.h>
36 #include <TSystem.h>
37 #include <TVirtualMC.h>
38
39 #include "AliMC.h"
40 #include "AliRun.h"
41 #include "AliMagF.h"
42 #include "AliConst.h"
43 #include "AliITSGeant3Geometry.h"
44 #include "AliTrackReference.h"
45 #include "AliITShit.h"
46 #include "AliITS.h"
47 #include "AliITSvSDD03.h"
48 #include "AliITSgeom.h"
49 #include "AliITSgeomSPD.h"
50 #include "AliITSgeomSDD.h"
51 #include "AliITSgeomSSD.h"
52 #include "AliITSDetType.h"
53 #include "AliITSresponseSPD.h"
54 #include "AliITSresponseSDD.h"
55 #include "AliITSresponseSSD.h"
56 #include "AliITSsegmentationSPD.h"
57 #include "AliITSsegmentationSDD.h"
58 #include "AliITSsegmentationSSD.h"
59 #include "AliITSsimulationSPD.h"
60 #include "AliITSsimulationSDD.h"
61 #include "AliITSsimulationSSD.h"
62 #include "AliITSClusterFinderSPD.h"
63 #include "AliITSClusterFinderSDD.h"
64 #include "AliITSClusterFinderSSD.h"
65
66
67 ClassImp(AliITSvSDD03)
68
69 //______________________________________________________________________
70 AliITSvSDD03::AliITSvSDD03() {
71     ////////////////////////////////////////////////////////////////////////
72     // Standard default constructor for the ITS SDD test beam 2002 version 1.
73     // Inputs:
74     //    none.
75     // Outputs:
76     //    none.
77     // Return:
78     //    A default created class.
79     ////////////////////////////////////////////////////////////////////////
80     Int_t i;
81
82     fIdN          = 0;
83     fIdName       = 0;
84     fIdSens       = 0;
85     fEuclidOut    = kFALSE; // Don't write Euclide file
86     fGeomDetOut   = kFALSE; // Don't write .det file
87     fGeomDetIn    = kFALSE; // Don't Read .det file
88     fMajorVersion = IsVersion();
89     fMinorVersion = -1;
90     for(i=0;i<60;i++) fRead[i] = '\0';
91     for(i=0;i<60;i++) fWrite[i] = '\0';
92     for(i=0;i<60;i++) fEuclidGeomDet[i] = '\0';
93 }
94 //______________________________________________________________________
95 AliITSvSDD03::AliITSvSDD03(const char *title) : AliITS("ITS", title){
96     ////////////////////////////////////////////////////////////////////////
97     //    Standard constructor for the ITS SDD testbeam 2002 version 1.
98     // Inputs:
99     //    const char *title    title for this ITS geometry.
100     // Outputs:
101     //    none.
102     // Return:
103     //    A standard created class.
104     ////////////////////////////////////////////////////////////////////////
105     Int_t i;
106
107     fIdN = 2;
108     fIdName = new TString[fIdN];
109     fIdName[0] = "IMBS";
110     fIdName[1] = "ITST";
111     fIdSens    = new Int_t[fIdN];
112     for(i=0;i<fIdN;i++) fIdSens[i] = 0;
113     fMajorVersion = IsVersion();
114     fMinorVersion = 2;
115     fEuclidOut    = kFALSE; // Don't write Euclide file
116     fGeomDetOut   = kFALSE; // Don't write .det file
117     fGeomDetIn    = kFALSE; // Don't Read .det file
118     SetThicknessDet1();
119     SetThicknessDet2();
120     SetThicknessChip1();
121     SetThicknessChip2();                         
122
123     fEuclidGeometry="$ALICE_ROOT/ITS/ITSgeometry_vSDD032.euc";
124     strncpy(fEuclidGeomDet,"$ALICE_ROOT/ITS/ITSgeometry_vSDD032.det",60);
125     strncpy(fRead,fEuclidGeomDet,60);
126     strncpy(fWrite,fEuclidGeomDet,60);
127 }
128 //______________________________________________________________________
129 AliITSvSDD03::AliITSvSDD03(const AliITSvSDD03 &source) : AliITS(source){
130     ////////////////////////////////////////////////////////////////////////
131     //     Copy Constructor for ITS SDD test beam 2002 version 1.
132     // This class is not to be copied. Function only dummy.
133     // Inputs:
134     //    const AliITSvSDD03 &source   The class to be copied
135     // Outputs:
136     //    none.
137     // Return:
138     //    A warning message.
139     ////////////////////////////////////////////////////////////////////////
140     if(&source == this) return;
141     Warning("Copy Constructor","Not allowed to copy AliITSvSDD03");
142     return;
143 }
144 //______________________________________________________________________
145 AliITSvSDD03& AliITSvSDD03::operator=(const AliITSvSDD03 &source){
146     ////////////////////////////////////////////////////////////////////////
147     //    Assignment operator for the ITS SDD test beam 2002 version 1.
148     // This class is not to be copied. Function only dummy.
149     // Inputs:
150     //    const AliITSvSDD03 &source   The class to be copied
151     // Outputs:
152     //    none.
153     // Return:
154     //    A Warning message
155     ////////////////////////////////////////////////////////////////////////
156     if(&source == this) return *this;
157     Warning("= operator","Not allowed to copy AliITSvSDD03");
158     return *this;
159 }
160 //______________________________________________________________________
161 AliITSvSDD03::~AliITSvSDD03() {
162     ////////////////////////////////////////////////////////////////////////
163     //    Standard destructor for the ITS SDD test beam 2002 version 1.
164     // Inputs:
165     //    none.
166     // Outputs:
167     //    none.
168     // Return:
169     //    none.
170     ////////////////////////////////////////////////////////////////////////
171 }
172 //______________________________________________________________________
173 void AliITSvSDD03::BuildGeometry(){
174     ////////////////////////////////////////////////////////////////////////
175     //    Geometry builder for the ITS SDD test beam 2002 version 1.
176     //    ALIC    ALICE Mother Volume
177     //     |- ITSV     ITS Mother Volume
178     //         |- IDET       Detector under Test
179     //         |   |- ITS0       SDD Si Chip
180     //         |   |  |- ITST      SDD Sensitivve Volume
181     //         |   |- IPC0 *5    Readout chip
182     //         |- ITEL *4    SDD Telescope
183     //             |- IMB0       SDD Si Chip
184     //             |   |- IMBS     SDD Sensitive volume
185     //             |- ICMB       Chip MiniBus.
186     // Inputs:
187     //    none.
188     // Outputs:
189     //    none.
190     // Return:
191     //    none.
192     ////////////////////////////////////////////////////////////////////////
193     // Get the top alice volume.
194     TNode *ALIC = gAlice->GetGeometry()->GetNode("alice");
195     ALIC->cd();
196
197     // Define ITS Mother Volume
198     Float_t data[3];
199     Float_t ddettest=200.0E-4,ddettelescope=300.0E-4;
200     Float_t dchipMiniBus=750.0E-4,dchiptest=300.0E-4;
201     //Float_t yposition= 0.0;
202     TRotMatrix *r0 = new TRotMatrix("ITSidrotm0","ITSidrotm0",
203                                     90.0,0,0.0,0,90.0,270.0);
204     data[0] = 10.0;
205     data[1] = 50.0;
206     data[2] = 100.0;
207     TBRIK *ITSVshape = new TBRIK("ITSVshape","ITS Logical Mother Volume","Air",
208                                  data[0],data[1],data[2]);
209     TNode *ITSV = new TNode("ITSV","ITS Mother Volume",ITSVshape,
210                             0.0,0.0,0.0,0,0);
211     ITSV->cd(); // set ourselve into ITSV subvolume of ALIC
212
213     // SDD part of telescope (MiniBuS)
214     data[0] = 0.705;
215     data[1] = 0.5*ddettelescope;
216     data[2] = 3.536;
217     TBRIK *IMB0shape = new TBRIK("IMB0shape","SDD wafer","Si",
218                                  data[0],data[1],data[2]);
219     Float_t detMiniBusX,detMiniBusY,detMiniBusZ;
220     data[0] = detMiniBusX = 0.64;
221     data[1] = detMiniBusY = 0.5*ddettelescope;
222     data[2] = detMiniBusZ = 3.48;
223     TBRIK *IMBSshape = new TBRIK("IMBSshape","SDD Sensitive volume","Si",
224                                  data[0],data[1],data[2]);
225     Float_t chipMiniBusX,chipMiniBusY,chipMiniBusZ;
226     data[0] = chipMiniBusX = 0.793;
227     data[1] = chipMiniBusY = 0.5*dchipMiniBus;
228     data[2] = chipMiniBusZ = 0.68;
229     TBRIK *ICMBshape = new TBRIK("ICMBshape","chip Minibus","Si",
230                                  data[0],data[1],data[2]);
231     data[0] = TMath::Max(detMiniBusX,chipMiniBusX);
232     data[1] = detMiniBusY+chipMiniBusY;
233     data[2] = TMath::Max(detMiniBusZ,chipMiniBusZ);
234     TBRIK *ITELshape = new TBRIK("ITELshape","ITELshape","Air",
235                                  data[0],data[1],data[2]);
236
237     // SDD under test
238     Float_t spdX,spdY,spdZ,spdchipX,spdchipY,spdchipZ;
239     data[0] = 0.705;
240     data[1] = ddettest;
241     data[2] = 3.536;
242     TBRIK *ITS0shape = new TBRIK("ITS0shape","SDD wafer","Si",
243                                  data[0],data[1],data[2]); // contains detector
244     data[0] = spdX = 0.64;
245     data[1] = spdY = ddettest;
246     data[2] = spdZ = 3.48;
247     TBRIK *ITSTshape = new TBRIK("ITSTshape","SDD sensitive volume","Si",
248                                  data[0],data[1],data[2]);
249     // ITS0 with no translation and unit rotation matrix.
250     data[0] = spdchipX = 0.793;
251     data[1] = spdchipY = dchiptest;
252     data[2] = spdchipZ = 0.68;
253     TBRIK *IPC0shape = new TBRIK("IPC0shape","Readout Chips","Si",
254                                  data[0],data[1],data[2]); // chip under test
255     data[0] = TMath::Max(spdchipX,spdX);
256     data[1] = spdY+spdchipY;
257     data[2] = TMath::Max(spdchipZ,spdZ);
258     TBRIK *IDETshape = new TBRIK("IDETshape","Detector Under Test","Air",
259                                  data[0],data[1],data[2]);
260     // Place volumes in geometry
261     Int_t i,j;
262     char name[20],title[50];
263     Double_t px=0.0,py=0.0,pz[4]={-38.0,0.0,0.0,0.0};
264     pz[1] = pz[0]+2.0;
265     pz[2] = pz[1]+38.0+spdY+spdchipY+34.5;
266     pz[3] = pz[2]+2.0;
267     TNode *ITEL[4],*ICMB[4],*IMB0[4],*IMBS[4];
268     TNode *IDET = new TNode("IDET","Detector Under Test",IDETshape,
269                             0.0,0.0,pz[1]+38.0,r0,0);
270     IDET->cd();
271     TNode *ITS0 = new TNode("ITS0","SDD Chip",ITS0shape,
272                             0.0,IDETshape->GetDy()-spdY,0.0,0,0);
273     TNode *IPC0[5];
274     for(i=0;i<5;i++) { //place readout chips on the back of SDD chip under test
275         sprintf(name,"IPC0%d",i);
276         sprintf(title,"Readout chip #%d",i+1);
277         j = i-2;
278         IPC0[i] = new TNode(name,title,IPC0shape,
279                             0.0,spdchipY-IDETshape->GetDy(),
280                             j*2.0*spdchipZ+j*0.25*(spdZ-5.*spdchipZ),0,0);
281     } // end for i
282     ITS0->cd();
283     TNode *ITST = new TNode("ITST","SDD sensitive volume",ITSTshape,
284                             0.0,0.0,0.0,0,0);
285     for(Int_t i=0;i<4;i++){
286         ITSV->cd();
287         sprintf(name,"ITEL%d",i);
288         sprintf(title,"Test beam telescope element #%d",i+1);
289         ITEL[i] = new TNode(name,title,ITELshape,px,py,pz[i],r0,0);
290         ITEL[i]->cd();
291         ICMB[i] = new TNode("ICMB","Chip MiniBus",ICMBshape,
292                             0.0,-ITELshape->GetDy()+detMiniBusY,0.0,0,0);
293         IMB0[i] = new TNode("IMB0","Chip MiniBus",IMB0shape,
294                             0.0, ITELshape->GetDy()-detMiniBusY,0.0,0,0);
295         IMB0[i]->cd();
296         IMBS[i] = new TNode("IMBS","IMBS",IMBSshape,0.0,0.0,0.0,0,0);
297         // place IMBS inside IMB0 with no translation and unit rotation matrix.
298     } // end for i
299     ALIC->cd();
300     ITST->SetLineColor(kYellow);
301     fNodes->Add(ITST);
302     for(i=0;i<4;i++){
303         IMBS[i]->SetLineColor(kGreen);
304         fNodes->Add(IMBS[i]);
305     } // end for i
306 }
307 //______________________________________________________________________
308 void AliITSvSDD03::CreateGeometry(){
309     ////////////////////////////////////////////////////////////////////////
310     //  This routine defines and Creates the geometry for version 1 of the ITS.
311     //    ALIC    ALICE Mother Volume
312     //     |- ITSV     ITS Mother Volume
313     //         |- IDET       Detector under Test (box containing SDD)
314     //         |   |-IDAI        Air inside box
315     //         |       |- ITS0       SDD Si Chip
316     //         |          |- ITST      SDD Sensitivve Volume
317     //         |- ITEL *10   SSD Telescope (plastic box containting SSD's)
318     //         |   |- ITAI       Air inside box
319     //         |       |- IMB0       SDD Si Chip
320     //         |           |- IMBS     SDD Sensitive volume
321     //         |-ISNT*4    Sintilator triggers
322     // Inputs:
323     //    none.
324     // Outputs:
325     //    none.
326     // Return:
327     //    none.
328     ////////////////////////////////////////////////////////////////////////
329     Float_t data[49];
330     // Define media off-set
331     Int_t *idtmed = fIdtmed->GetArray()+1; // array of media indexes
332     Int_t idrotm[4]; // Array of rotation matrix indexes
333     //Float_t ddettest=200.0E-4,ddettelescope=300.0E-4;
334     //Float_t dchipMiniBus=750.0E-4,dchiptest=300.0E-4;
335     //Float_t yposition= 0.0;
336     const Float_t kmm=0.1,kcm=1.0,kmicm=0.001;
337
338     // Define Rotation-reflextion Matrixes needed
339     // 0 is the unit matrix
340     AliMatrix(idrotm[0], 90.0,0.0, 0.0,0.0, 90.0,270.0); // SDD and SSD X
341     AliMatrix(idrotm[1], 90.0,0.0, 0.0,0.0, 90.0,270.0); // SSD Y
342     data[0] = 100.0*kmm;
343     data[1] = 100.0*kmm;
344     data[2] = 800.0*kmm;
345     gMC->Gsvolu("ITSV","BOX ",idtmed[0],data,3);
346     gMC->Gspos("ITSV",1,"ALIC",0.0,0.0,0.0,0,"ONLY");
347
348     //cout << "idtmed[0]=" << idtmed[0]<<endl;
349     //cout << "idtmed[1]=" << idtmed[1]<<endl;
350     // Crossed sintilator triggers (2 in front 2 in back)
351     AliMatrix(idrotm[2],90.0,0.0,90.0,90.0,90.0,0.0);//Rotate about Z 90 degree
352     data[0] = 10.0*kcm;
353     data[1] = 2.0*kcm;
354     data[2] = 2.0*kmm;
355     gMC->Gsvolu("ISNT","BOX ",idtmed[3],data,3);
356     gMC->Gspos("ISNT",1,"ITSV",0.0,0.0,800.0*kmm+data[2],0,"ONLY");
357     gMC->Gspos("ISNT",2,"ITSV",0.0,0.0,800.0*kmm,idrotm[2],"ONLY");
358     gMC->Gspos("ISNT",3,"ITSV",0.0,0.0,-800.0*kmm,0,"ONLY");
359     gMC->Gspos("ISNT",4,"ITSV",0.0,0.0,-800.0*kmm-data[2],idrotm[2],"ONLY");
360     Float_t detMiniBusX,detMiniBusY,detMiniBusZ;
361     // SSD part of telescope (MiniBuS)
362     data[0] = detMiniBusX = 10600.0*kmicm;
363     data[1] = detMiniBusY = 0.0150*kcm;
364     data[2] = detMiniBusZ = 1.1*kcm;
365     gMC->Gsvolu("IMB0", "BOX ", idtmed[1], data, 3);   // contains detector
366     data[0] = 0.5*384*50*kmicm;
367     data[1] = 0.1499*kcm;
368     data[2] = 1.0*kcm;
369     gMC->Gsvolu("IMBS","BOX ",idtmed[1],data,3); // sensitive detecor volulme
370     gMC->Gspos("IMBS",1,"IMB0",0.0,0.0,0.0,0,"ONLY"); // place IMBS inside
371     // Box containing SSD's
372     data[0] = 4.0*kcm;
373     data[1] = 0.500*kcm;
374     data[2] = 5.0*kcm;
375     gMC->Gsvolu("ITAI","BOX ",idtmed[0],data,3);
376     // Plastic box size = insize + thickness.
377     data[0] = data[0] + 2.0*kmm;
378     data[1] = data[1] + 200.0*kmicm;
379     data[2] = data[2] + 2.0*kmm;
380     gMC->Gsvolu("ITEL","BOX ",idtmed[4],data,3);
381     gMC->Gspos("ITAI",1,"ITEL",0.0,0.0,0.0,0,"ONLY");
382     gMC->Gspos("IMB0",1,"ITAI",0.0,data[1]-detMiniBusY,0.0,0,"ONLY");
383
384     // SDD under test
385     Float_t sddX,sddY,sddZ;
386     data[0] = sddX = 3.62500*kcm;
387     data[1] = sddY = 0.01500*kcm;
388     data[2] = sddZ = 4.37940*kcm;
389     gMC->Gsvolu("ITS0", "BOX ", idtmed[1], data, 3);   // contains detector
390     data[0] = 3.50860*kcm;
391     data[1] = 0.01499*kcm;
392     data[2] = 3.76320*kcm;
393     gMC->Gsvolu("ITST","BOX ",idtmed[1],data,3);// sensitive detecor volume
394     gMC->Gspos("ITST",1,"ITS0",0.0,0.0,0.0,0,"ONLY"); // place ITST inside
395     // Box containing SDD under test
396     data[0] = 4.0*kcm;
397     data[1] = 0.5*kcm;
398     data[2] = 5.0*kcm;
399     gMC->Gsvolu("IDAI","BOX ",idtmed[0],data,3);
400     data[0] = data[0] + 2.0*kmm;
401     data[1] = data[1] + 200.0*kmicm;
402     data[2] = data[2] = 2.0*kmm;
403     gMC->Gsvolu("IDET","BOX ",idtmed[0],data,3);
404     gMC->Gspos("IDAI",1,"IDET",0.0,0.0,0.0,0,"ONLY");
405     gMC->Gspos("ITS0",1,"IDAI",0.0,0.0,0.0,0,"ONLY");
406
407     // Positions detectors, Beam Axis Z, X to the right, Y up to the sky.
408     Float_t p00X,p00Y,p00Z,p01X,p01Y,p01Z,p10X,p10Y,p10Z,p11X,p11Y,p11Z;
409     p00X = 0.0*kcm;
410     p00Y = 0.0*kcm;
411     p00Z = -694*kmm;
412     gMC->Gspos("ITEL",1,"ITSV",p00X,p00Y,p00Z,idrotm[0],"ONLY");//SSD X
413     p01X = 0.0*kcm;
414     p01Y = 0.0*kcm;
415     p01Z = -684*kmm;
416     gMC->Gspos("ITEL",2,"ITSV",p01X,p01Y,p01Z,idrotm[1],"ONLY");//SSD Y
417     p01X = 0.0*kcm;
418     p01Y = 0.0*kcm;
419     p01Z = -612*kmm;
420     gMC->Gspos("ITEL",3,"ITSV",p01X,p01Y,p01Z,idrotm[0],"ONLY");//SSD X
421     p01X = 0.0*kcm;
422     p01Y = 0.0*kcm;
423     p01Z = -602*kmm;
424     Float_t pdetX,pdetY,pdetZ;
425     gMC->Gspos("ITEL",4,"ITSV",p01X,p01Y,p01Z,idrotm[1],"ONLY");//SSD Y
426     pdetX = 0.0*kcm;
427     pdetY = 0.0*kcm;
428     pdetZ = 0.0*kcm;
429     gMC->Gspos("IDET",1,"ITSV",pdetX,pdetY,pdetZ,idrotm[0],"ONLY");// Detecor
430     p10X = 0.0*kcm;
431     p10Y = 0.0*kcm;
432     p10Z = +450.0*kmm;
433     gMC->Gspos("ITEL",5,"ITSV",p10X,p10Y,p10Z,idrotm[0],"ONLY");//SSD X
434     p11X = 0.0*kcm;
435     p11Y = 0.0*kcm;
436     p11Z = +460.0*kcm;
437     gMC->Gspos("ITEL",6,"ITSV",p11X,p11Y,p11Z,idrotm[1],"ONLY");//SSD Y
438     p11X = 0.0*kcm;
439     p11Y = 0.0*kcm;
440     p11Z = +540.0*kcm;
441     gMC->Gspos("ITEL",7,"ITSV",p11X,p11Y,p11Z,idrotm[0],"ONLY");//SSD X
442     p11X = 0.0*kcm;
443     p11Y = 0.0*kcm;
444     p11Z = +550.0*kcm;
445     gMC->Gspos("ITEL",8,"ITSV",p11X,p11Y,p11Z,idrotm[1],"ONLY");//SSD Y
446     p11X = 0.0*kcm;
447     p11Y = 0.0*kcm;
448     p11Z = +737.0*kcm;
449     gMC->Gspos("ITEL",9,"ITSV",p11X,p11Y,p11Z,idrotm[0],"ONLY");//SSD X
450     p11X = 0.0*kcm;
451     p11Y = 0.0*kcm;
452     p11Z = +747.0*kcm;
453     gMC->Gspos("ITEL",10,"ITSV",p11X,p11Y,p11Z,idrotm[1],"ONLY");//SSD Y
454 }
455 //______________________________________________________________________
456 void AliITSvSDD03::CreateMaterials(){
457     ////////////////////////////////////////////////////////////////////////
458     //
459     // Create ITS SDD test beam materials
460     //     This function defines the default materials used in the Geant
461     // Monte Carlo simulations for the geometries AliITSv1, AliITSv3,
462     // AliITSvSDD03.
463     // In general it is automatically replaced by
464     // the CreatMaterials routine defined in AliITSv?. Should the function
465     // CreateMaterials not exist for the geometry version you are using this
466     // one is used. See the definition found in AliITSv5 or the other routine
467     // for a complete definition.
468     //
469     // Inputs:
470     //    none.
471     // Outputs:
472     //    none.
473     // Return:
474     //    none.
475     /////////////////////////////////////////////////////////////////////////
476     Float_t tmaxfdSi = 0.1; // Degree
477     Float_t stemaxSi = 0.0075; // cm
478     Float_t deemaxSi = 0.1; // Fraction of particle's energy 0<deemax<=1
479     Float_t epsilSi  = 1.0E-4;//
480     Float_t stminSi  = 0.0; // cm "Default value used"
481
482     Float_t tmaxfdAir = 0.1; // Degree
483     Float_t stemaxAir = .10000E+01; // cm
484     Float_t deemaxAir = 0.1; // Fraction of particle's energy 0<deemax<=1
485     Float_t epsilAir  = 1.0E-4;//
486     Float_t stminAir  = 0.0; // cm "Default value used"
487     Int_t   ifield = gAlice->Field()->Integ();
488     Float_t fieldm = gAlice->Field()->Max();
489     //
490     const Float_t kgpcm3=1.0,kcm=1.0;
491     //
492     Float_t z[10],a[10],w[10];
493
494     z[0] = 7.0; a[0] = 14.00674; w[0] = 0.80;
495     z[1] = 8.0; a[1] = 15.99940; w[1] = 0.20;
496     AliMixture(1,"AIR$",a,z,0.12050E-02*kgpcm3,2,w);
497     AliMedium(1,"AIR$",1,0,ifield,fieldm,tmaxfdAir,stemaxAir,deemaxAir,
498               epsilAir,stminAir);
499
500     AliMaterial(2,"SI$",0.28086E+02,0.14000E+02,0.23300E+01,
501                 0.93600E+01*kcm,0.99900E+03);
502     AliMedium(2,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
503               epsilSi,stminSi);
504     // sintilator is basicaly polystyrene
505     z[0] = 7.0; a[0] = 14.00674; w[0] = 0.80;
506     AliMixture(3,"Sintilator$",a,z,1.032*kgpcm3,0,w);
507     AliMedium(3,"Sintilator$",3,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
508               epsilSi,stminSi);
509     // assumed to be Lucite/Plexiglas
510     z[0] = 7.0; a[0] = 14.00674; w[0] = 0.80;
511     AliMixture(4,"PlasticBox$",a,z,1.18*kgpcm3,0,w);
512     AliMedium(4,"PlasticBox$",4,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
513               epsilSi,stminSi);
514 }
515 //______________________________________________________________________
516 void AliITSvSDD03::InitAliITSgeom(){
517     //     Based on the geometry tree defined in Geant 3.21, this
518     // routine initilizes the Class AliITSgeom from the Geant 3.21 ITS geometry
519     // sturture.
520     // Inputs:
521     //    none.
522     // Outputs:
523     //    none.
524     // Return:
525     //    none.
526     if(strcmp(gMC->GetName(),"TGeant3")) {
527         Error("InitAliITSgeom",
528                 "Wrong Monte Carlo. InitAliITSgeom uses TGeant3 calls");
529         return;
530     } // end if
531     cout << "Reading Geometry transformation directly from Geant 3." << endl;
532     const Int_t np=384;
533     const Float_t pitch=50.E-4;/*cm*/
534     Float_t box[3]={0.5*pitch*(Float_t)np,150.E-4,1.0},p[np],n[np];
535     const Int_t ltypess = 2;
536     const Int_t nlayers = 6;
537     const Int_t ndeep = 6;
538     Int_t itsGeomTreeNames[ltypess][ndeep],lnam[20],lnum[20];
539     Int_t nlad[nlayers],ndet[nlayers];
540     Double_t t[3],r[10];
541     Float_t  par[20],att[20];
542     Int_t    npar,natt,idshape,imat,imed;
543     AliITSGeant3Geometry *ig = new AliITSGeant3Geometry();
544     Int_t mod,typ,lay,lad,det,cpy,i,j,k;
545     Char_t names[ltypess][ndeep][4];
546     Int_t itsGeomTreeCopys[ltypess][ndeep];
547     Char_t *namesA[ltypess][ndeep] = {
548      {"ALIC","ITSV","ITEL","ITIA","IMB0","IMBS"}, // lay=5
549      {"ALIC","ITSV","IDET","IDAI","ITS0","ITST"}};// Test SDD
550     Int_t itsGeomTreeCopysA[ltypess][ndeep]= {{1,1,10,1,1,1},// lay=5
551                                               {1,1,1,1,1,1}};//lay=3 TestSDD
552     for(i=0;i<ltypess;i++)for(j=0;j<ndeep;j++){
553         for(k=0;k<4;k++) names[i][j][k] = namesA[i][j][k];
554         itsGeomTreeCopys[i][j] = itsGeomTreeCopysA[i][j];
555     } // end for i,j
556     for(i=0;i<np;i++){// Fill in anode and cathode strip locations (lower edge)
557         p[i] = 0.5*pitch*(Float_t)np + pitch*(Float_t)i;
558         n[i] = pitch*(Float_t)np - p[i];
559     } // end for i
560     // Sorry, but this is not very pritty code. It should be replaced
561     // at some point with a version that can search through the geometry
562     // tree its self.
563     cout << "Reading Geometry informaton from Geant3 common blocks" << endl;
564     for(i=0;i<20;i++) lnam[i] = lnum[i] = 0;
565     for(i=0;i<ltypess;i++)for(j=0;j<ndeep;j++) 
566         strncpy((char*) &itsGeomTreeNames[i][j],names[i][j],4);
567     //  itsGeomTreeNames[i][j] = ig->StringToInt(names[i][j]);
568     mod = 5;
569     if(fITSgeom!=0) delete fITSgeom;
570     nlad[0]=1;nlad[1]=1;nlad[2]=1;nlad[3]=1;nlad[4]=1;
571     ndet[0]=1;ndet[1]=1;ndet[2]=1;ndet[3]=1;ndet[4]=1;
572     fITSgeom = new AliITSgeom(0,nlayers,nlad,ndet,mod);
573     for(typ=1;typ<=ltypess;typ++){
574         for(j=0;j<ndeep;j++) lnam[j] = itsGeomTreeNames[typ-1][j];
575         for(j=0;j<ndeep;j++) lnum[j] = itsGeomTreeCopys[typ-1][j];
576         lad = 1;
577         det = 1;
578         for(cpy=1;cpy<=itsGeomTreeCopys[typ-1][2];cpy++){
579             lnum[2] = cpy;
580             lay = cpy;
581             if(cpy>2 && typ==1) lay = cpy +1;
582             if(typ==2) lay = 3;
583             mod = lay-1;
584             ig->GetGeometry(ndeep,lnam,lnum,t,r,idshape,npar,natt,par,att,
585                             imat,imed);
586             if(!(fITSgeom->IsShapeDefined((Int_t)kSDD)))if(typ==1){
587                 fITSgeom->CreatMatrix(mod,lay,lad,det,kSDD,t,r);
588                 fITSgeom->ReSetShape(kSDD,new AliITSgeomSDD256(npar,par));
589             }else{
590                 fITSgeom->CreatMatrix(mod,lay,lad,det,kSSD,t,r);
591                 fITSgeom->ReSetShape(kSSD,new AliITSgeomSSD(box,0.0,0.0,
592                                                             np+1,p,np+1,n));
593             } // end if
594         } // end for cpy
595     } // end for typ
596     return;
597 }
598 //______________________________________________________________________
599 void AliITSvSDD03::Init(){
600     ////////////////////////////////////////////////////////////////////////
601     //     Initialise the ITS after it has been created.
602     // Inputs:
603     //    none.
604     // Outputs:
605     //    none.
606     // Return:
607     //    none.
608     ////////////////////////////////////////////////////////////////////////
609     Int_t i;
610
611     cout << endl;
612     for(i=0;i<26;i++) cout << "*";
613     cout << " AliITSvSDD03" << fMinorVersion << "_Init ";
614     for(i=0;i<25;i++) cout << "*";cout << endl;
615 //
616     if(fRead[0]=='\0') strncpy(fRead,fEuclidGeomDet,60);
617     if(fWrite[0]=='\0') strncpy(fWrite,fEuclidGeomDet,60);
618     if(fITSgeom!=0) delete fITSgeom;
619     fITSgeom = new AliITSgeom();
620     if(fGeomDetIn) fITSgeom->ReadNewFile(fRead);
621     if(!fGeomDetIn) this->InitAliITSgeom();
622     if(fGeomDetOut) fITSgeom->WriteNewFile(fWrite);
623     AliITS::Init();
624 //
625     for(i=0;i<72;i++) cout << "*";
626     cout << endl;
627     fIDMother = gMC->VolId("ITSV"); // ITS Mother Volume ID.
628 }
629 //______________________________________________________________________
630 void AliITSvSDD03::SetDefaults(){
631     // sets the default segmentation, response, digit and raw cluster classes
632     // Inputs:
633     //    none.
634     // Outputs:
635     //    none.
636     // Return:
637     //    none.
638     const Float_t kconv = 1.0e+04; // convert cm to microns
639
640     Info("SetDefaults","Setting up only SDD detector");
641
642     AliITSDetType *iDetType;
643     AliITSgeomSDD *s1;
644     AliITSgeomSSD *s2;
645     //Int_t i;
646     //Float_t bx[256],bz[280];
647
648     // SDD
649     iDetType=DetType(kSDD);
650     s1 = (AliITSgeomSDD*) fITSgeom->GetShape(kSDD);// Get shape info. Do it this way for now.
651     AliITSresponseSDD *resp1=new AliITSresponseSDD("simulated");
652     SetResponseModel(kSDD,resp1);
653     AliITSsegmentationSDD *seg1=new AliITSsegmentationSDD(fITSgeom,resp1);
654     seg1->SetDetSize(s1->GetDx()*kconv, // base this on AliITSgeomSDD
655                      s1->GetDz()*2.*kconv, // for now.
656                      s1->GetDy()*2.*kconv); // x,z,y full width in microns.
657     seg1->SetNPads(256,256);// Use AliITSgeomSDD for now
658     SetSegmentationModel(kSDD,seg1);
659     const char *kData1=(iDetType->GetResponseModel())->DataType();
660     const char *kopt=iDetType->GetResponseModel()->ZeroSuppOption();
661     if((!strstr(kopt,"2D")) && (!strstr(kopt,"1D")) || strstr(kData1,"real") ){
662         iDetType->ClassNames("AliITSdigit","AliITSRawClusterSDD");
663     } else iDetType->ClassNames("AliITSdigitSDD","AliITSRawClusterSDD");
664 //    SetSimulationModel(kSDD,new AliITSsimulationSDD(seg1,resp1));
665 //    iDetType->ReconstructionModel(new AliITSClusterFinderSDD());
666
667     // SSD  Layer 5
668     iDetType=DetType(kSSD);
669     s2 = (AliITSgeomSSD*) fITSgeom->GetShape(kSSD);// Get shape info. Do it this way for now.
670     AliITSresponse *resp2=new AliITSresponseSSD("simulated");
671     SetResponseModel(kSSD,resp2);
672     AliITSsegmentationSSD *seg2=new AliITSsegmentationSSD(fITSgeom);
673     seg2->SetDetSize(s2->GetDx()*2.*kconv, // base this on AliITSgeomSSD
674                      s2->GetDz()*2.*kconv, // for now.
675                      s2->GetDy()*2.*kconv); // x,z,y full width in microns.
676     seg2->SetPadSize(95.,0.); // strip x pitch in microns
677     seg2->SetNPads(768,0); // number of strips on each side.
678     seg2->SetAngles(0.0075,0.0275); // strip angels rad P and N side.
679     seg2->SetAnglesLay5(0.0075,0.0275); // strip angels rad P and N side.
680     seg2->SetAnglesLay6(0.0275,0.0075); // strip angels rad P and N side.
681     SetSegmentationModel(kSSD,seg2); 
682     const char *kData2=(iDetType->GetResponseModel())->DataType();
683     if(strstr(kData2,"real") ) iDetType->ClassNames("AliITSdigit",
684                                                     "AliITSRawClusterSSD");
685     else iDetType->ClassNames("AliITSdigitSSD","AliITSRawClusterSSD");
686 //    SetSimulationModel(kSSD,new AliITSsimulationSSD(seg2,resp2));
687 //    iDetType->ReconstructionModel(new AliITSClusterFinderSSD());
688
689     if(kNTYPES>3){
690         Warning("SetDefaults",
691                 "Only the four basic detector types are initialised!");
692     }// end if
693     return;
694 }
695 //______________________________________________________________________
696 void AliITSvSDD03::DrawModule(){
697     ////////////////////////////////////////////////////////////////////////
698     //     Draw a shaded view of the ITS SDD test beam version 1.
699     // Inputs:
700     //    none.
701     // Outputs:
702     //    none.
703     // Return:
704     //    none.
705     ////////////////////////////////////////////////////////////////////////
706     // Set everything unseen
707     gMC->Gsatt("*", "seen", -1);
708     // Set ALIC mother visible
709     gMC->Gsatt("ALIC","SEEN",0);
710     // Set ALIC ITS visible
711     gMC->Gsatt("ITSV","SEEN",0);
712     // Set ALIC Telescopes visible
713     gMC->Gsatt("ITEL","SEEN",0);
714     // Set ALIC detetcor visible
715     gMC->Gsatt("IDET","SEEN",0);
716     // Set Detector chip mother visible and drawn
717     gMC->Gsatt("IPC0","SEEN",1);
718     // Set Detector mother visible and drawn
719     gMC->Gsatt("ITS0","SEEN",1);
720     // Set minibus chip mother visible and drawn
721     gMC->Gsatt("ICMB","SEEN",1);
722     // Set minibus mother visible and drawn
723     gMC->Gsatt("IMB0","SEEN",1);
724 }
725 //______________________________________________________________________
726 void AliITSvSDD03::StepManager(){
727     ////////////////////////////////////////////////////////////////////////
728     //    Called for every step in the ITS SDD test beam, then calles the 
729     // AliITShit class  creator with the information to be recoreded about
730     //  that hit.
731     //     The value of the macro ALIITSPRINTGEOM if set to 1 will allow the
732     // printing of information to a file which can be used to create a .det
733     // file read in by the routine CreateGeometry(). If set to 0 or any other
734     // value except 1, the default behavior, then no such file is created nor
735     // it the extra variables and the like used in the printing allocated.
736     // Inputs:
737     //    none.
738     // Outputs:
739     //    none.
740     // Return:
741     //    none.
742     ////////////////////////////////////////////////////////////////////////
743     Int_t         copy, id;
744     TLorentzVector position, momentum;
745     static TLorentzVector position0;
746     static Int_t stat0=0;
747     if((id=gMC->CurrentVolID(copy) == fIDMother)&&
748        (gMC->IsTrackEntering()||gMC->IsTrackExiting())){
749         copy = fTrackReferences->GetEntriesFast();
750         TClonesArray &lTR = *fTrackReferences;
751         // Fill TrackReference structure with this new TrackReference.
752         new(lTR[copy]) AliTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber());
753     } // if Outer ITS mother Volume
754     if(!(this->IsActive())){
755         return;
756     } // end if !Active volume.
757     Int_t   vol[5];
758     TClonesArray &lhits = *fHits;
759     //
760     // Track status
761     vol[3] = 0;
762     vol[4] = 0;
763     if(gMC->IsTrackInside())      vol[3] +=  1;
764     if(gMC->IsTrackEntering())    vol[3] +=  2;
765     if(gMC->IsTrackExiting())     vol[3] +=  4;
766     if(gMC->IsTrackOut())         vol[3] +=  8;
767     if(gMC->IsTrackDisappeared()) vol[3] += 16;
768     if(gMC->IsTrackStop())        vol[3] += 32;
769     if(gMC->IsTrackAlive())       vol[3] += 64;
770     //
771     // Fill hit structure.
772     if(!(gMC->TrackCharge())) return;
773     id = gMC->CurrentVolID(copy);
774     if(id==fIdSens[0]){  // Volume name "IMBS"
775         vol[2] = 1; // detector
776         id = gMC->CurrentVolOffID(3,copy);
777         //detector copy in the ladder = 1<->4  (ITS1 < I101 < I103 < I10A)
778         vol[1] = copy; // ladder
779         vol[0] = 5; // Lay
780         if(copy>4){
781             vol[0] = 6;
782             vol[1] = copy - 4;
783         } // end if
784         if(copy>2) vol[0]++;
785     } else if(id == fIdSens[1]){ // Volume name "ITST"
786         vol[0] = 3; // layer
787         vol[1] = 1; // ladder
788         id = gMC->CurrentVolOffID(3,copy);
789         //detector copy in the ladder = 1<->4  (ITS2 < I1D1 < I1D3 < I20A)
790         vol[2] = 1;  // detector
791     } else return; // end if
792     //
793     gMC->TrackPosition(position);
794     gMC->TrackMomentum(momentum);
795     vol[4] = stat0;
796     if(gMC->IsTrackEntering()){
797         position0 = position;
798         stat0 = vol[3];
799         return;
800     } // end if IsEntering
801     // Fill hit structure with this new hit only for non-entrerance hits.
802     else new(lhits[fNhits++]) AliITShit(fIshunt,
803                                         gAlice->GetMCApp()->GetCurrentTrackNumber(),vol,
804                                         gMC->Edep(),gMC->TrackTime(),position,
805                                         position0,momentum);
806     //
807     position0 = position;
808     stat0 = vol[3];
809
810     return;
811 }
812