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