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