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