]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSvSPD02.cxx
Correct material budget for tracks stopped at TRD (S.Arcelli)
[u/mrichter/AliRoot.git] / ITS / AliITSvSPD02.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 /* $Id$ */
17
18 #include <Riostream.h>
19 #include <TMath.h>
20 #include <TGeometry.h>
21 #include <TNode.h>
22 #include <TLorentzVector.h>
23 #include <TClonesArray.h>
24 #include <TBRIK.h>
25
26 #include "AliRun.h"
27 #include "AliMagF.h"
28 #include "AliITSGeant3Geometry.h"
29 #include "AliTrackReference.h"
30 #include "AliITShit.h"
31 #include "AliITS.h"
32 #include "AliITSvSPD02.h"
33 #include "AliITSgeom.h"
34 #include "AliITSgeomSPD.h"
35 #include "AliITSgeomSSD.h"
36 #include "AliITSDetType.h"
37 #include "AliITSresponseSPD.h"
38 #include "AliITSresponseSDD.h"
39 #include "AliITSresponseSSD.h"
40 #include "AliITSsegmentationSPD.h"
41 #include "AliITSsegmentationSDD.h"
42 #include "AliITSsegmentationSSD.h"
43 #include "AliITSsimulationSPDdubna.h"
44 #include "AliITSsimulationSDD.h"
45 #include "AliITSsimulationSSD.h"
46 #include "AliMC.h"
47
48 ///////////////////////////////////////////////////////////////////////
49 // Step manager and 
50 // geometry class
51 // for the ITS 
52 // SPD test beam
53 // geometry of summer 2002
54 // 
55 ///////////////////////////////////////////////////////////////////////
56 ClassImp(AliITSvSPD02)
57
58 //______________________________________________________________________
59 AliITSvSPD02::AliITSvSPD02():
60 AliITS(),
61 fGeomDetOut(kFALSE),
62 fGeomDetIn(kFALSE),
63 fMajorVersion(1),
64 fMinorVersion(2),
65 fGeomNumber(2002),
66 fEuclidGeomDet(),
67 fRead(),
68 fWrite(),
69 fDet1(300.0),
70 fDet2(300.0),
71 fChip1(300.0),
72 fChip2(300.0){
73     ////////////////////////////////////////////////////////////////////////
74     // Standard default constructor for the ITS SPD 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     for(i=0;i<60;i++) fRead[i] = '\0';
85     for(i=0;i<60;i++) fWrite[i] = '\0';
86     for(i=0;i<60;i++) fEuclidGeomDet[i] = '\0';
87 }
88 //______________________________________________________________________
89 AliITSvSPD02::AliITSvSPD02(const char *title,Int_t gn) : AliITS("ITS", title),
90 fGeomDetOut(kFALSE),
91 fGeomDetIn(kFALSE),
92 fMajorVersion(1),
93 fMinorVersion(2),
94 fGeomNumber(2002),
95 fEuclidGeomDet(),
96 fRead(),
97 fWrite(),
98 fDet1(300.0),
99 fDet2(300.0),
100 fChip1(300.0),
101 fChip2(300.0){
102     ////////////////////////////////////////////////////////////////////////
103     //    Standard constructor for the ITS SPD testbeam 2002 version 1.
104     // Inputs:
105     //    const char *title    title for this ITS geometry.
106     //    Int_t      gn        Geometry version number (year) default 2002.
107     // Outputs:
108     //    none.
109     // Return:
110     //    A standard created class.
111     ////////////////////////////////////////////////////////////////////////
112     Int_t i;
113
114     fGeomNumber = gn;
115     fIdN = 2;
116     fIdName = new TString[fIdN];
117     fIdName[0] = "IMBS";
118     fIdName[1] = "ITST";
119     fIdSens    = new Int_t[fIdN];
120     for(i=0;i<fIdN;i++) fIdSens[i] = 0;
121     SetThicknessDet1();
122     SetThicknessDet2();
123     SetThicknessChip1();
124     SetThicknessChip2();                         
125
126     fEuclidGeometry="$ALICE_ROOT/ITS/ITSgeometry_vSPD022.euc";
127     strncpy(fEuclidGeomDet,"$ALICE_ROOT/ITS/ITSgeometry_vSPD022.det",60);
128     strncpy(fRead,fEuclidGeomDet,60);
129     strncpy(fWrite,fEuclidGeomDet,60);
130 }
131 //______________________________________________________________________
132 AliITSvSPD02::AliITSvSPD02(const AliITSvSPD02 &source) : AliITS(source){
133     ////////////////////////////////////////////////////////////////////////
134     //     Copy Constructor for ITS SPD test beam 2002 version 1.
135     // This class is not to be copied. Function only dummy.
136     // Inputs:
137     //    const AliITSvSPD02 &source   The class to be copied
138     // Outputs:
139     //    none.
140     // Return:
141     //    A warning message.
142     ////////////////////////////////////////////////////////////////////////
143     if(&source == this) return;
144     Warning("Copy Constructor","Not allowed to copy AliITSvSPD02");
145     return;
146 }
147 //______________________________________________________________________
148 AliITSvSPD02& AliITSvSPD02::operator=(const AliITSvSPD02 &source){
149     ////////////////////////////////////////////////////////////////////////
150     //    Assignment operator for the ITS SPD test beam 2002 version 1.
151     // This class is not to be copied. Function only dummy.
152     // Inputs:
153     //    const AliITSvSPD02 &source   The class to be copied
154     // Outputs:
155     //    none.
156     // Return:
157     //    A Warning message
158     ////////////////////////////////////////////////////////////////////////
159     if(&source == this) return *this;
160     Warning("= operator","Not allowed to copy AliITSvSPD02");
161     return *this;
162 }
163 //______________________________________________________________________
164 AliITSvSPD02::~AliITSvSPD02() {
165     ////////////////////////////////////////////////////////////////////////
166     //    Standard destructor for the ITS SPD test beam 2002 version 1.
167     // Inputs:
168     //    none.
169     // Outputs:
170     //    none.
171     // Return:
172     //    none.
173     ////////////////////////////////////////////////////////////////////////
174 }
175 //______________________________________________________________________
176 void AliITSvSPD02::BuildGeometry(){
177     ////////////////////////////////////////////////////////////////////////
178     //    Geometry builder for the ITS SPD test beam 2002 version 1.
179     //    ALIC    ALICE Mother Volume
180     //     |- ITSV     ITS Mother Volume
181     //         |- IDET       Detector under Test
182     //         |   |- ITS0       SPD Si Chip
183     //         |   |  |- ITST      SPD Sensitivve Volume
184     //         |   |- IPC0 *5    Readout chip
185     //         |- ITEL *4    SPD Telescope
186     //             |- IMB0       SPD Si Chip
187     //             |   |- IMBS     SPD Sensitive volume
188     //             |- ICMB       Chip MiniBus.
189     // Inputs:
190     //    none.
191     // Outputs:
192     //    none.
193     // Return:
194     //    none.
195     ////////////////////////////////////////////////////////////////////////
196     // Get the top alice volume.
197
198     switch (fGeomNumber){
199     case 2002:
200         BuildGeometry2002();
201         break;
202     default:
203         BuildGeometry2002();
204         break;
205     } // end switch
206 }
207 //______________________________________________________________________
208 void AliITSvSPD02::BuildGeometry2002(){
209     ////////////////////////////////////////////////////////////////////////
210     //    Geometry builder for the ITS SPD test beam 2002 version 1.
211     //    ALIC    ALICE Mother Volume
212     //     |- ITSV     ITS Mother Volume
213     //         |- IDET       Detector under Test
214     //         |   |- ITS0       SPD Si Chip
215     //         |   |  |- ITST      SPD Sensitivve Volume
216     //         |   |- IPC0 *5    Readout chip
217     //         |- ITEL *4    SPD Telescope
218     //             |- IMB0       SPD Si Chip
219     //             |   |- IMBS     SPD Sensitive volume
220     //             |- ICMB       Chip MiniBus.
221     // Inputs:
222     //    none.
223     // Outputs:
224     //    none.
225     // Return:
226     //    none.
227     ////////////////////////////////////////////////////////////////////////
228     // Get the top alice volume.
229     TNode *aALIC = gAlice->GetGeometry()->GetNode("alice");
230     aALIC->cd();
231
232     // Define ITS Mother Volume
233     Float_t data[3];
234     Float_t ddettest=200.0E-4,ddettelescope=300.0E-4;
235     Float_t dchipMiniBus=750.0E-4,dchiptest=300.0E-4;
236     //Float_t yposition= 0.0;
237     TRotMatrix *r0 = new TRotMatrix("ITSidrotm0","ITSidrotm0",
238                                     90.0,0,0.0,0,90.0,270.0);
239     data[0] = 10.0;
240     data[1] = 50.0;
241     data[2] = 100.0;
242     TBRIK *iITSVshape = new TBRIK("ITSVshape",
243                                   "ITS Logical Mother Volume","Air",
244                                   data[0],data[1],data[2]);
245     TNode *iITSV = new TNode("ITSV","ITS Mother Volume",iITSVshape,
246                              0.0,0.0,0.0,0,0);
247     iITSV->cd(); // set ourselve into ITSV subvolume of aALIC
248
249     // SPD part of telescope (MiniBuS)
250     data[0] = 0.705;
251     data[1] = 0.5*ddettelescope;
252     data[2] = 3.536;
253     TBRIK *iIMB0shape = new TBRIK("IMB0shape","SPD wafer","Si",
254                                  data[0],data[1],data[2]);
255     Float_t detMiniBusX,detMiniBusY,detMiniBusZ;
256     data[0] = detMiniBusX = 0.64;
257     data[1] = detMiniBusY = 0.5*ddettelescope;
258     data[2] = detMiniBusZ = 3.48;
259     TBRIK *iIMBSshape = new TBRIK("IMBSshape","SPD Sensitive volume","Si",
260                                  data[0],data[1],data[2]);
261     Float_t chipMiniBusX,chipMiniBusY,chipMiniBusZ;
262     data[0] = chipMiniBusX = 0.793;
263     data[1] = chipMiniBusY = 0.5*dchipMiniBus;
264     data[2] = chipMiniBusZ = 0.68;
265     TBRIK *iICMBshape = new TBRIK("ICMBshape","chip Minibus","Si",
266                                  data[0],data[1],data[2]);
267     data[0] = TMath::Max(detMiniBusX,chipMiniBusX);
268     data[1] = detMiniBusY+chipMiniBusY;
269     data[2] = TMath::Max(detMiniBusZ,chipMiniBusZ);
270     TBRIK *iITELshape = new TBRIK("ITELshape","ITELshape","Air",
271                                  data[0],data[1],data[2]);
272
273     // SPD under test
274     Float_t spdX,spdY,spdZ,spdchipX,spdchipY,spdchipZ;
275     data[0] = 0.705;
276     data[1] = ddettest;
277     data[2] = 3.536;
278     TBRIK *iITS0shape = new TBRIK("ITS0shape","SPD wafer","Si",
279                                  data[0],data[1],data[2]); // contains detector
280     data[0] = spdX = 0.64;
281     data[1] = spdY = ddettest;
282     data[2] = spdZ = 3.48;
283     TBRIK *iITSTshape = new TBRIK("ITSTshape","SPD sensitive volume","Si",
284                                  data[0],data[1],data[2]);
285     // ITS0 with no translation and unit rotation matrix.
286     data[0] = spdchipX = 0.793;
287     data[1] = spdchipY = dchiptest;
288     data[2] = spdchipZ = 0.68;
289     TBRIK *iIPC0shape = new TBRIK("IPC0shape","Readout Chips","Si",
290                                  data[0],data[1],data[2]); // chip under test
291     data[0] = TMath::Max(spdchipX,spdX);
292     data[1] = spdY+spdchipY;
293     data[2] = TMath::Max(spdchipZ,spdZ);
294     TBRIK *iIDETshape = new TBRIK("IDETshape","Detector Under Test","Air",
295                                  data[0],data[1],data[2]);
296     // Place volumes in geometry
297     Int_t i,j;
298     char name[20],title[50];
299     Double_t px=0.0,py=0.0,pz[4]={-38.0,0.0,0.0,0.0};
300     pz[1] = pz[0]+2.0;
301     pz[2] = pz[1]+38.0+spdY+spdchipY+34.5;
302     pz[3] = pz[2]+2.0;
303     TNode *iITEL[4],*iICMB[4],*iIMB0[4],*iIMBS[4];
304     TNode *iIDET = new TNode("IDET","Detector Under Test",iIDETshape,
305                             0.0,0.0,pz[1]+38.0,r0,0);
306     iIDET->cd();
307     TNode *iITS0 = new TNode("ITS0","SPD Chip",iITS0shape,
308                             0.0,iIDETshape->GetDy()-spdY,0.0,0,0);
309     TNode *iIPC0[5];
310     for(i=0;i<5;i++) { //place readout chips on the back of SPD chip under test
311         sprintf(name,"IPC0%d",i);
312         sprintf(title,"Readout chip #%d",i+1);
313         j = i-2;
314         iIPC0[i] = new TNode(name,title,iIPC0shape,
315                              0.0,spdchipY-iIDETshape->GetDy(),
316                              j*2.0*spdchipZ+j*0.25*(spdZ-5.*spdchipZ),0,0);
317     } // end for i
318     iITS0->cd();
319     TNode *iITST = new TNode("ITST","SPD sensitive volume",iITSTshape,
320                              0.0,0.0,0.0,0,0);
321     for(Int_t i=0;i<4;i++){
322         iITSV->cd();
323         sprintf(name,"ITEL%d",i);
324         sprintf(title,"Test beam telescope element #%d",i+1);
325         iITEL[i] = new TNode(name,title,iITELshape,px,py,pz[i],r0,0);
326         iITEL[i]->cd();
327         iICMB[i] = new TNode("ICMB","Chip MiniBus",iICMBshape,
328                              0.0,-iITELshape->GetDy()+detMiniBusY,0.0,0,0);
329         iIMB0[i] = new TNode("IMB0","Chip MiniBus",iIMB0shape,
330                              0.0, iITELshape->GetDy()-detMiniBusY,0.0,0,0);
331         iIMB0[i]->cd();
332         iIMBS[i] = new TNode("IMBS","IMBS",iIMBSshape,0.0,0.0,0.0,0,0);
333         // place IMBS inside IMB0 with no translation and unit rotation matrix.
334     } // end for i
335     aALIC->cd();
336     iITST->SetLineColor(kYellow);
337     fNodes->Add(iITST);
338     for(i=0;i<4;i++){
339         iIMBS[i]->SetLineColor(kGreen);
340         fNodes->Add(iIMBS[i]);
341     } // end for i
342 }
343 //______________________________________________________________________
344 void AliITSvSPD02::CreateGeometry(){
345     ////////////////////////////////////////////////////////////////////////
346     //  This routine defines and Creates the geometry for version 1 of the ITS.
347     //    ALIC    ALICE Mother Volume
348     //     |- ITSV     ITS Mother Volume
349     //         |- IDET       Detector under Test
350     //         |   |- ITS0       SPD Si Chip
351     //         |   |  |- ITST      SPD Sensitivve Volume
352     //         |   |- IPC0 *5    Readout chip
353     //         |- ITEL *4    SPD Telescope
354     //             |- IMB0       SPD Si Chip
355     //             |   |- IMBS     SPD Sensitive volume
356     //             |- ICMB       Chip MiniBus.
357     // Inputs:
358     //    none.
359     // Outputs:
360     //    none.
361     // Return:
362     //    none.
363     ////////////////////////////////////////////////////////////////////////
364
365     switch (fGeomNumber){
366     case 2002:
367         CreateGeometry2002();
368         break;
369     default:
370         CreateGeometry2002();
371         break;
372     } // end switch
373 }
374 //______________________________________________________________________
375 void AliITSvSPD02::CreateGeometry2002(){
376     ////////////////////////////////////////////////////////////////////////
377     //  This routine defines and Creates the geometry for version 1 of the ITS.
378     //    ALIC    ALICE Mother Volume
379     //     |- ITSV     ITS Mother Volume
380     //         |- IDET       Detector under Test
381     //         |   |- ITS0       SPD Si Chip
382     //         |   |  |- ITST      SPD Sensitivve Volume
383     //         |   |- IPC0 *5    Readout chip
384     //         |- ITEL *4    SPD Telescope
385     //             |- IMB0       SPD Si Chip
386     //             |   |- IMBS     SPD Sensitive volume
387     //             |- ICMB       Chip MiniBus.
388     // Inputs:
389     //    none.
390     // Outputs:
391     //    none.
392     // Return:
393     //    none.
394     ////////////////////////////////////////////////////////////////////////
395     Float_t data[49];
396     // Define media off-set
397     Int_t *idtmed = fIdtmed->GetArray()+1; // array of media indexes
398     Int_t idrotm[4]; // Array of rotation matrix indexes
399     Float_t ddettest=200.0E-4,ddettelescope=300.0E-4;
400     Float_t dchipMiniBus=750.0E-4,dchiptest=300.0E-4;
401     Float_t yposition= 0.0;
402
403     if(gMC==0) return;
404     // Define Rotation-reflextion Matrixes needed
405     // 0 is the unit matrix
406     AliMatrix(idrotm[0], 90.0,0.0, 0.0,0.0, 90.0,270.0);
407     data[0] = 10.0;
408     data[1] = 50.0;
409     data[2] = 100.0;
410     gMC->Gsvolu("ITSV","BOX",idtmed[0],data,3);
411     gMC->Gspos("ITSV",1,"ALIC",0.0,0.0,0.0,0,"ONLY");
412
413     //cout << "idtmed[0]=" << idtmed[0]<<endl;
414     //cout << "idtmed[1]=" << idtmed[1]<<endl;
415     Float_t detMiniBusX,detMiniBusY,detMiniBusZ;
416     // SPD part of telescope (MiniBuS)
417     data[0] = detMiniBusX = 0.705;
418     data[1] = detMiniBusY = 0.5*ddettelescope;
419     data[2] = detMiniBusZ = 3.536;
420     gMC->Gsvolu("IMB0", "BOX ", idtmed[1], data, 3);   // contains detector
421     data[0] = 0.64;
422     data[1] = 0.5*ddettelescope;
423     data[2] = 3.48;
424     gMC->Gsvolu("IMBS","BOX ",idtmed[1],data,3); // sensitive detecor volulme
425     gMC->Gspos("IMBS",1,"IMB0",0.0,0.0,0.0,0,"ONLY"); // place IMBS inside
426     // IMB0 with no translation and unit rotation matrix.
427     Float_t chipMiniBusX,chipMiniBusY,chipMiniBusZ;
428     data[0] = chipMiniBusX = 0.793;
429     data[1] = chipMiniBusY = 0.5*dchipMiniBus;
430     data[2] = chipMiniBusZ = 0.68;
431     gMC->Gsvolu("ICMB","BOX ",idtmed[1],data, 3);   // chip Minibus
432     data[0] = TMath::Max(detMiniBusX,chipMiniBusX);
433     data[1] = detMiniBusY+chipMiniBusY;
434     data[2] = TMath::Max(detMiniBusZ,chipMiniBusZ);
435     gMC->Gsvolu("ITEL","BOX ",idtmed[0],data,3);
436     gMC->Gspos("IMB0",1,"ITEL",0.0,data[1]-detMiniBusY,0.0,0,"ONLY");
437     gMC->Gspos("ICMB",1,"ITEL",0.0,-data[1]+chipMiniBusY,0.0,0,"ONLY");
438
439     // SPD under test
440     Float_t spdX,spdY,spdZ,spdchipX,spdchipY,spdchipZ;
441     data[0] = spdX = 0.705;
442     data[1] = spdY = 0.5*ddettest;
443     data[2] = spdZ = 3.536;
444     gMC->Gsvolu("ITS0", "BOX ", idtmed[1], data, 3);   // contains detector
445     data[0] = 0.64;
446     data[1] = 0.5*ddettest;
447     data[2] = 3.48;
448     gMC->Gsvolu("ITST","BOX ",idtmed[1],data,3);// sensitive detecor volume
449     gMC->Gspos("ITST",1,"ITS0",0.0,0.0,0.0,0,"ONLY"); // place ITST inside
450     // ITS0 with no translation and unit rotation matrix.
451     data[0] = spdchipX = 0.793;
452     data[1] = spdchipY = 0.5*dchiptest;
453     data[2] = spdchipZ = 0.68;
454     gMC->Gsvolu("IPC0", "BOX ", idtmed[1],data,3);   // chip under test
455     data[0] = TMath::Max(spdchipX,spdX);
456     data[1] = spdY+spdchipY;
457     data[2] = TMath::Max(spdchipZ,spdZ);
458     gMC->Gsvolu("IDET","BOX ",idtmed[0],data,3);
459     gMC->Gspos("ITS0",1,"IDET",0.0,data[1]-spdY,0.0,0,"ONLY");
460     for(Int_t i=-2;i<3;i++) gMC->Gspos("IPC0",i+3,"IDET",0.0,-data[1]+spdchipY,
461               i*2.*spdchipZ+i*0.25*(spdZ-5.*spdchipZ),0,"ONLY");
462
463     // Positions detectors, Beam Axis Z, X to the right, Y up to the sky.
464     Float_t p00X,p00Y,p00Z,p01X,p01Y,p01Z,p10X,p10Y,p10Z,p11X,p11Y,p11Z;
465     p00X = 0.0;
466     p00Y = 0.0;
467     p00Z = -38.0;
468     gMC->Gspos("ITEL",1,"ITSV",p00X,p00Y,p00Z,idrotm[0],"ONLY");
469     p01X = 0.0;
470     p01Y = 0.0;
471     p01Z = p00Z+2.0;
472     gMC->Gspos("ITEL",2,"ITSV",p01X,p01Y,p01Z,idrotm[0],"ONLY");
473     Float_t pdetX,pdetY,pdetZ;
474     pdetX = 0.0;
475     pdetY = 0.0+yposition;
476     pdetZ = p01Z+38.0;
477     gMC->Gspos("IDET",1,"ITSV",pdetX,pdetY,pdetZ,idrotm[0],"ONLY");
478     p10X = 0.0;
479     p10Y = 0.0;
480     p10Z = pdetZ + 34.5;
481     gMC->Gspos("ITEL",3,"ITSV",p10X,p10Y,p10Z,idrotm[0],"ONLY");
482     p11X = 0.0;
483     p11Y = 0.0;
484     p11Z = p10Z+2.0;
485     gMC->Gspos("ITEL",4,"ITSV",p11X,p11Y,p11Z,idrotm[0],"ONLY");
486 }
487 //______________________________________________________________________
488 void AliITSvSPD02::CreateMaterials(){
489     ////////////////////////////////////////////////////////////////////////
490     //
491     // Create ITS SPD test beam materials
492     //     This function defines the default materials used in the Geant
493     // Monte Carlo simulations for the geometries AliITSv1, AliITSv3,
494     // AliITSvSPD02.
495     // In general it is automatically replaced by
496     // the CreatMaterials routine defined in AliITSv?. Should the function
497     // CreateMaterials not exist for the geometry version you are using this
498     // one is used. See the definition found in AliITSv5 or the other routine
499     // for a complete definition.
500     //
501     // Inputs:
502     //    none.
503     // Outputs:
504     //    none.
505     // Return:
506     //    none.
507     /////////////////////////////////////////////////////////////////////////
508
509     switch (fGeomNumber){
510     case 2002:
511         CreateMaterials2002();
512         break;
513     default:
514         CreateMaterials2002();
515         break;
516     } // end switch
517 }
518 //______________________________________________________________________
519 void AliITSvSPD02::CreateMaterials2002(){
520     ////////////////////////////////////////////////////////////////////////
521     //
522     // Create ITS SPD test beam materials
523     //     This function defines the default materials used in the Geant
524     // Monte Carlo simulations for the geometries AliITSv1, AliITSv3,
525     // AliITSvSPD02.
526     // In general it is automatically replaced by
527     // the CreatMaterials routine defined in AliITSv?. Should the function
528     // CreateMaterials not exist for the geometry version you are using this
529     // one is used. See the definition found in AliITSv5 or the other routine
530     // for a complete definition.
531     //
532     // Inputs:
533     //    none.
534     // Outputs:
535     //    none.
536     // Return:
537     //    none.
538     /////////////////////////////////////////////////////////////////////////
539     Float_t tmaxfdSi = 0.1; // Degree
540     Float_t stemaxSi = 0.0075; // cm
541     Float_t deemaxSi = 0.1; // Fraction of particle's energy 0<deemax<=1
542     Float_t epsilSi  = 1.0E-4;//
543     Float_t stminSi  = 0.0; // cm "Default value used"
544
545     Float_t tmaxfdAir = 0.1; // Degree
546     Float_t stemaxAir = .10000E+01; // cm
547     Float_t deemaxAir = 0.1; // Fraction of particle's energy 0<deemax<=1
548     Float_t epsilAir  = 1.0E-4;//
549     Float_t stminAir  = 0.0; // cm "Default value used"
550     Int_t   ifield = gAlice->Field()->Integ();
551     Float_t fieldm = gAlice->Field()->Max();
552
553     AliMaterial(1,"AIR$",0.14610E+02,0.73000E+01,0.12050E-02,
554                 0.30423E+05,0.99900E+03);
555     AliMedium(1,"AIR$",1,0,ifield,fieldm,tmaxfdAir,stemaxAir,deemaxAir,
556               epsilAir,stminAir);
557
558     AliMaterial(2,"SI$",0.28086E+02,0.14000E+02,0.23300E+01,
559                 0.93600E+01,0.99900E+03);
560     AliMedium(2,"SI$",2,0,ifield,fieldm,tmaxfdSi,stemaxSi,deemaxSi,
561               epsilSi,stminSi);
562 }
563 //______________________________________________________________________
564 void AliITSvSPD02::InitAliITSgeom(){
565     //     Based on the geometry tree defined in Geant 3.21, this
566     // routine initilizes the Class AliITSgeom from the Geant 3.21 ITS geometry
567     // sturture.
568     // Inputs:
569     //    none.
570     // Outputs:
571     //    none.
572     // Return:
573     //    none.
574     const Int_t kltypess = 2;
575     const Int_t knlayers = 5;
576     const Int_t kndeep = 5;
577     Int_t itsGeomTreeNames[kltypess][kndeep],lnam[20],lnum[20];
578     Int_t nlad[knlayers],ndet[knlayers];
579     Double_t t[3],r[10];
580     Float_t  par[20],att[20];
581     Int_t    npar,natt,idshape,imat,imed;
582     AliITSGeant3Geometry *ig=0;
583     Int_t mod=0,typ=0,lay=0,lad=0,det=0,cpy=0,i=0,j=0,k=0;
584
585     if(gMC==0) {// No Monti Carlo to init. Default set fITSgeom by hand
586         if(fITSgeom!=0) delete fITSgeom;
587         nlad[0]=1;nlad[1]=1;nlad[2]=1;nlad[3]=1;nlad[4]=1;
588         ndet[0]=1;ndet[1]=1;ndet[2]=1;ndet[3]=1;ndet[4]=1;
589         fITSgeom = new AliITSgeom(0,knlayers,nlad,ndet,mod);
590         r[0] = 1.0; r[1] =  0.0; r[2] = 0.0;
591         r[3] = 0.0; r[4] =  0.0; r[5] = 1.0;
592         r[6] = 0.0; r[7] = -1.0; r[8] = 0.0; r[9] = 1.0; // not Unit.
593         Double_t tt[5][3]={{0.0,0.0,-37.9625}, // 0,0,P00Z+.5*chipMiniBusY
594                            {0.0,0.0,-35.9625}, // 0,0,P00Z+2+.5*chipMiniBusY
595                            {0.0,0.0,  2.0150}, // 0,0,P01Z+38+SPDchipY
596                            {0.0,0.0, 36.5375}, //0,0,PdetZ+34.5+.5*chipMiniBusY
597                            {0.0,0.0, 38.5375}};// 0,0,P00Z+2+.5*chipMiniBusY
598         for(mod=0;mod<5;mod++){
599             lay = 1;
600             lad = 1;
601             det = mod+1;
602             t[0] = tt[mod][0]; t[1] = tt[mod][1]; t[2] = tt[mod][2];
603             fITSgeom->CreatMatrix(mod,lay,lad,det,kSPD,t,r);
604             npar=3;par[0]=0.64;par[1]=0.5*300.0E-4;par[2]=3.48;
605             fITSgeom->ReSetShape(kSPD,new AliITSgeomSPD425Short(npar,par));
606         } // end for det
607         return;
608     } // end if gMC==0
609     if(strcmp(gMC->GetName(),"TGeant3")) {
610         Error("InitAliITSgeom",
611               "Wrong Monte Carlo. InitAliITSgeom uses TGeant3 calls");
612         return;
613     } // end if
614     cout << "Reading Geometry transformation directly from Geant 3." << endl;
615     ig = new AliITSGeant3Geometry();
616     Char_t names[kltypess][kndeep][4];
617     Int_t itsGeomTreeCopys[kltypess][kndeep];
618     const char *namesA[kltypess][kndeep] = {
619         {"ALIC","ITSV","ITEL","IMB0","IMBS"}, // lay=1
620         {"ALIC","ITSV","IDET","ITS0","ITST"}};// Test SPD
621     Int_t itsGeomTreeCopysA[kltypess][kndeep]= {{1,1,4,1,1},// lay=1
622                                                 {1,1,1,1,1}};//lay=2 TestSPD
623     for(i=0;i<kltypess;i++)for(j=0;j<kndeep;j++){
624         for(k=0;k<4;k++) names[i][j][k] = namesA[i][j][k];
625         itsGeomTreeCopys[i][j] = itsGeomTreeCopysA[i][j];
626     } // end for i,j
627     // Sorry, but this is not very pritty code. It should be replaced
628     // at some point with a version that can search through the geometry
629     // tree its self.
630     cout << "Reading Geometry informaton from Geant3 common blocks" << endl;
631     for(i=0;i<20;i++) lnam[i] = lnum[i] = 0;
632     for(i=0;i<kltypess;i++)for(j=0;j<kndeep;j++) 
633         strncpy((char*) &itsGeomTreeNames[i][j],names[i][j],4);
634     //  itsGeomTreeNames[i][j] = ig->StringToInt(names[i][j]);
635     mod = 5;
636     if(fITSgeom!=0) delete fITSgeom;
637     nlad[0]=1;nlad[1]=1;nlad[2]=1;nlad[3]=1;nlad[4]=1;
638     ndet[0]=1;ndet[1]=1;ndet[2]=1;ndet[3]=1;ndet[4]=1;
639     fITSgeom = new AliITSgeom(0,knlayers,nlad,ndet,mod);
640     for(typ=1;typ<=kltypess;typ++){
641         for(j=0;j<kndeep;j++) lnam[j] = itsGeomTreeNames[typ-1][j];
642         for(j=0;j<kndeep;j++) lnum[j] = itsGeomTreeCopys[typ-1][j];
643         lad = 1;
644         det = 1;
645         for(cpy=1;cpy<=itsGeomTreeCopys[typ-1][2];cpy++){
646             lnum[2] = cpy;
647             lay = cpy;
648             if(cpy>2 && typ==1) lay = cpy +1;
649             if(typ==2) lay = 3;
650             mod = lay-1;
651             ig->GetGeometry(kndeep,lnam,lnum,t,r,idshape,npar,natt,par,att,
652                             imat,imed);
653             fITSgeom->CreatMatrix(mod,lay,lad,det,kSPD,t,r);
654             if(!(fITSgeom->IsShapeDefined((Int_t)kSPD)))
655                 fITSgeom->ReSetShape(kSPD,
656                                      new AliITSgeomSPD425Short(npar,par));
657         } // end for cpy
658     } // end for typ
659     return;
660 }
661 //______________________________________________________________________
662 void AliITSvSPD02::Init(){
663     ////////////////////////////////////////////////////////////////////////
664     //     Initialise the ITS after it has been created.
665     // Inputs:
666     //    none.
667     // Outputs:
668     //    none.
669     // Return:
670     //    none.
671     ////////////////////////////////////////////////////////////////////////
672     Int_t i;
673
674     cout << endl;
675     for(i=0;i<26;i++) cout << "*";
676     cout << " ITSvSPD02" << fMinorVersion << "_Init ";
677     for(i=0;i<25;i++) cout << "*";cout << endl;
678
679     if(fRead[0]=='\0') strncpy(fRead,fEuclidGeomDet,60);
680     if(fWrite[0]=='\0') strncpy(fWrite,fEuclidGeomDet,60);
681     if(fITSgeom!=0) delete fITSgeom;
682     fITSgeom = new AliITSgeom();
683     if(fGeomDetIn) fITSgeom->ReadNewFile(fRead);
684     if(!fGeomDetIn) this->InitAliITSgeom();
685     if(fGeomDetOut) fITSgeom->WriteNewFile(fWrite);
686     AliITS::Init();
687
688     for(i=0;i<72;i++) cout << "*";
689     cout << endl;
690     if(gMC) fIDMother = gMC->VolId("ITSV"); // ITS Mother Volume ID.
691     else fIDMother = 0;
692 }
693 //______________________________________________________________________
694 void AliITSvSPD02::SetDefaults(){
695     // sets the default segmentation, response, digit and raw cluster classes
696     // Inputs:
697     //    none.
698     // Outputs:
699     //    none.
700     // Return:
701     //    none.
702     const Float_t kconv = 1.0e+04; // convert cm to microns
703
704     Info("SetDefaults","Setting up only SPD detector");
705
706     AliITSDetType *iDetType;
707     AliITSgeomSPD  *s0;
708     Int_t i;
709     Float_t bx[256],bz[280];
710
711     //SPD
712     iDetType=DetType(kSPD);
713     // Get shape info. Do it this way for now.
714     s0 = (AliITSgeomSPD*) fITSgeom->GetShape(kSPD);
715     AliITSresponse *resp0=new AliITSresponseSPD();
716     resp0->SetTemperature();
717     resp0->SetDistanceOverVoltage();
718     SetResponseModel(kSPD,resp0);
719     AliITSsegmentationSPD *seg0=new AliITSsegmentationSPD(fITSgeom);
720     seg0->SetDetSize(s0->GetDx()*2.*kconv, // base this on AliITSgeomSPD
721                      s0->GetDz()*2.*kconv, // for now.
722                      s0->GetDy()*2.*kconv); // x,z,y full width in microns.
723     seg0->SetNPads(256,160);// Number of Bins in x and z
724     for(i=000;i<256;i++) bx[i] =  50.0; // in x all are 50 microns.
725     for(i=000;i<160;i++) bz[i] = 425.0; // most are 425 microns except below
726     for(i=160;i<280;i++) bz[i] =   0.0; // Outside of detector.
727     bz[ 31] = bz[ 32] = 625.0; // first chip boundry
728     bz[ 63] = bz[ 64] = 625.0; // first chip boundry
729     bz[ 95] = bz[ 96] = 625.0; // first chip boundry
730     bz[127] = bz[128] = 625.0; // first chip boundry
731     bz[160] = 425.0; // Set so that there is no zero pixel size for fNz.
732     seg0->SetBinSize(bx,bz); // Based on AliITSgeomSPD for now.
733     SetSegmentationModel(kSPD,seg0);
734     // set digit and raw cluster classes to be used
735     const char *kData0=(iDetType->GetResponseModel())->DataType();
736     if (strstr(kData0,"real")) iDetType->ClassNames("AliITSdigit",
737                                                     "AliITSRawClusterSPD");
738     else iDetType->ClassNames("AliITSdigitSPD","AliITSRawClusterSPD");
739 //    SetSimulationModel(kSPD,new AliITSsimulationSPDdubna(seg0,resp0));
740 //    iDetType->ReconstructionModel(new AliITSClusterFinderSPD());
741     delete DetType(kSDD);
742     fDetTypes->AddAt(0,kSDD);
743     delete DetType(kSSD);
744     fDetTypes->AddAt(0,kSSD);
745 /*
746     SetResponseModel(kSDD,new AliITSresponseSDD());
747     SetSegmentationModel(kSDD,new AliITSsegmentationSDD());
748     DetType(kSDD)->ClassNames("AliITSdigitSDD","AliITSRawClusterSDD");
749
750     SetResponseModel(kSSD,new AliITSresponseSSD());
751     SetSegmentationModel(kSSD,new AliITSsegmentationSSD());
752     DetType(kSSD)->ClassNames("AliITSdigitSSD","AliITSRawClusterSSD");
753 */
754     if(kNTYPES>3){
755         Warning("SetDefaults",
756                 "Only the four basic detector types are initialised!");
757     }// end if
758     return;
759 }
760 //______________________________________________________________________
761 void AliITSvSPD02::SetDefaultSimulation(){
762     // sets the default simulation.
763     // Inputs:
764     //      none.
765     // Outputs:
766     //      none.
767     // Return:
768     //      none.
769
770     AliITSDetType *iDetType;
771     AliITSsimulation *sim;
772     AliITSsegmentation *seg;
773     AliITSresponse *res;
774     iDetType=DetType(kSPD);
775     if(iDetType){
776         sim = iDetType->GetSimulationModel();
777         if (!sim) {
778             seg = (AliITSsegmentation*)iDetType->GetSegmentationModel();
779             res = (AliITSresponse*)iDetType->GetResponseModel();
780             sim = new AliITSsimulationSPDdubna(seg,res,0);
781             SetSimulationModel(kSPD,sim);
782         }else{ // simulation exists, make sure it is set up properly.
783             if(sim->GetResponseModel()==0) sim->SetResponseModel(
784                 (AliITSresponse*)iDetType->GetResponseModel());
785             if(sim->GetSegmentationModel()==0) sim->SetSegmentationModel(
786                 (AliITSsegmentation*)iDetType->GetSegmentationModel());
787             sim->Init();
788         } // end if
789     } // end if iDetType
790     iDetType=DetType(kSDD);
791     if(iDetType){
792         sim = iDetType->GetSimulationModel();
793         if (!sim) {
794             seg = (AliITSsegmentation*)iDetType->GetSegmentationModel();
795             res = (AliITSresponse*)iDetType->GetResponseModel();
796             sim = new AliITSsimulationSDD(seg,res);
797             SetSimulationModel(kSDD,sim);
798         }else{ // simulation exists, make sure it is set up properly.
799             if(sim->GetResponseModel()==0) sim->SetResponseModel(
800                 (AliITSresponse*)iDetType->GetResponseModel());
801             if(sim->GetSegmentationModel()==0) sim->SetSegmentationModel(
802                 (AliITSsegmentation*)iDetType->GetSegmentationModel());
803             sim->Init();
804         } //end if
805     } // end if iDetType
806     iDetType=DetType(kSSD);
807     if(iDetType){
808         sim = iDetType->GetSimulationModel();
809         if (!sim) {
810             seg = (AliITSsegmentation*)iDetType->GetSegmentationModel();
811             res = (AliITSresponse*)iDetType->GetResponseModel();
812             sim = new AliITSsimulationSSD(seg,res);
813             SetSimulationModel(kSSD,sim);
814         }else{ // simulation exists, make sure it is set up properly.
815             if(sim->GetResponseModel()==0) sim->SetResponseModel(
816                 (AliITSresponse*)iDetType->GetResponseModel());
817             if(sim->GetSegmentationModel()==0) sim->SetSegmentationModel(
818                 (AliITSsegmentation*)iDetType->GetSegmentationModel());
819             sim->Init();
820         } // end if
821     } // end if iDetType
822 }
823 //______________________________________________________________________
824 void AliITSvSPD02::DrawModule() const {
825     ////////////////////////////////////////////////////////////////////////
826     //     Draw a shaded view of the ITS SPD test beam version 1.
827     // Inputs:
828     //    none.
829     // Outputs:
830     //    none.
831     // Return:
832     //    none.
833     ////////////////////////////////////////////////////////////////////////
834     // Set everything unseen
835     gMC->Gsatt("*", "seen", -1);
836     // Set ALIC mother visible
837     gMC->Gsatt("ALIC","SEEN",0);
838     // Set ALIC ITS visible
839     gMC->Gsatt("ITSV","SEEN",0);
840     // Set ALIC Telescopes visible
841     gMC->Gsatt("ITEL","SEEN",0);
842     // Set ALIC detetcor visible
843     gMC->Gsatt("IDET","SEEN",0);
844     // Set Detector chip mother visible and drawn
845     gMC->Gsatt("IPC0","SEEN",1);
846     // Set Detector mother visible and drawn
847     gMC->Gsatt("ITS0","SEEN",1);
848     // Set minibus chip mother visible and drawn
849     gMC->Gsatt("ICMB","SEEN",1);
850     // Set minibus mother visible and drawn
851     gMC->Gsatt("IMB0","SEEN",1);
852 }
853 //______________________________________________________________________
854 void AliITSvSPD02::StepManager(){
855     ////////////////////////////////////////////////////////////////////////
856     //    Called for every step in the ITS SPD test beam, then calles the 
857     // AliITShit class  creator with the information to be recoreded about
858     //  that hit.
859     //     The value of the macro ALIITSPRINTGEOM if set to 1 will allow the
860     // printing of information to a file which can be used to create a .det
861     // file read in by the routine CreateGeometry(). If set to 0 or any other
862     // value except 1, the default behavior, then no such file is created nor
863     // it the extra variables and the like used in the printing allocated.
864     // Inputs:
865     //    none.
866     // Outputs:
867     //    none.
868     // Return:
869     //    none.
870     ////////////////////////////////////////////////////////////////////////
871     Int_t         copy, id;
872     TLorentzVector position, momentum;
873     static TLorentzVector position0;
874     static Int_t stat0=0;
875     if((id=gMC->CurrentVolID(copy) == fIDMother)&&
876        (gMC->IsTrackEntering()||gMC->IsTrackExiting())){
877         copy = fTrackReferences->GetEntriesFast();
878         TClonesArray &lTR = *fTrackReferences;
879         // Fill TrackReference structure with this new TrackReference.
880         new(lTR[copy]) AliTrackReference(gAlice->GetMCApp()->
881                                          GetCurrentTrackNumber());
882     } // if Outer ITS mother Volume
883     if(!(this->IsActive())){
884         return;
885     } // end if !Active volume.
886     Int_t   vol[5];
887     TClonesArray &lhits = *fHits;
888     //
889     // Track status
890     vol[3] = 0;
891     vol[4] = 0;
892     if(gMC->IsTrackInside())      vol[3] +=  1;
893     if(gMC->IsTrackEntering())    vol[3] +=  2;
894     if(gMC->IsTrackExiting())     vol[3] +=  4;
895     if(gMC->IsTrackOut())         vol[3] +=  8;
896     if(gMC->IsTrackDisappeared()) vol[3] += 16;
897     if(gMC->IsTrackStop())        vol[3] += 32;
898     if(gMC->IsTrackAlive())       vol[3] += 64;
899     //
900     // Fill hit structure.
901     if(!(gMC->TrackCharge())) return;
902     id = gMC->CurrentVolID(copy);
903     if(id==fIdSens[0]){  // Volume name "IMBS"
904         vol[2] = vol[1] = 1; // Det, ladder
905         id = gMC->CurrentVolOffID(2,copy);
906         //detector copy in the ladder = 1<->4  (ITS1 < I101 < I103 < I10A)
907         vol[0] = copy; // Lay
908         if(copy>2) vol[0]++;
909     } else if(id == fIdSens[1]){ // Volume name "ITST"
910         vol[0] = 3; // layer
911         vol[1] = 1; // ladder
912         id = gMC->CurrentVolOffID(2,copy);
913         //detector copy in the ladder = 1<->4  (ITS2 < I1D1 < I1D3 < I20A)
914         vol[2] = 1;  // detector
915     } else return; // end if
916     //
917     gMC->TrackPosition(position);
918     gMC->TrackMomentum(momentum);
919     vol[4] = stat0;
920     if(gMC->IsTrackEntering()){
921         position0 = position;
922         stat0 = vol[3];
923         return;
924     } // end if IsEntering
925     // Fill hit structure with this new hit only for non-entrerance hits.
926     else new(lhits[fNhits++]) AliITShit(fIshunt,
927                                   gAlice->GetMCApp()->GetCurrentTrackNumber(),
928                                         vol,gMC->Edep(),gMC->TrackTime(),
929                                         position,position0,momentum);
930     //
931     position0 = position;
932     stat0 = vol[3];
933
934     return;
935 }
936