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