]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSv5.cxx
Fixed bug in distructor. Can't distroy fixxed length arrays. Thanks Peter.
[u/mrichter/AliRoot.git] / ITS / AliITSv5.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.29  2001/02/09 00:05:31  nilsen
19 Added fMajor/MinorVersion variables and made other changes to better make
20 use of the new code changes in AliITSgeom related classes.
21
22 Revision 1.28  2001/02/02 23:57:28  nilsen
23 Added include file that are no londer included in AliITSgeom.h
24
25 Revision 1.27  2001/01/30 09:23:13  hristov
26 Streamers removed (R.Brun)
27
28 Revision 1.26  2000/11/30 11:13:11  barbera
29  Added changes suggested by Federico Carminati on nov, 30, 2000
30
31 Revision 1.25  2000/10/05 20:50:00  nilsen
32 Now using root generated streamers.
33
34 Revision 1.14.4.12  2000/10/02 16:04:03  barbera
35 Forward declarations added
36
37 Revision 1.22  2000/07/10 16:07:19  fca
38 Release version of ITS code
39
40 Revision 1.14.4.4  2000/05/19 10:10:21  nilsen
41 fix for bug with HP and Sun unix + fix for event display in ITS-working branch
42
43 Revision 1.14.4.3  2000/03/04 23:46:38  nilsen
44 Fixed up the comments/documentation.
45
46 Revision 1.14.4.2  2000/03/02 21:53:02  nilsen
47 To make it compatable with the changes in AliRun/AliModule.
48
49 Revision 1.14.4.1  2000/01/12 19:03:33  nilsen
50 This is the version of the files after the merging done in December 1999.
51 See the ReadMe110100.txt file for details
52
53 Revision 1.14  1999/10/22 08:16:49  fca
54 Correct destructors, thanks to I.Hrivnacova
55
56 Revision 1.13  1999/10/06 10:15:19  fca
57 Correct bug in allocation of layer name and add destructor
58
59 Revision 1.12  1999/10/05 08:05:09  fca
60 Minor corrections for uninitialised variables.
61
62 Revision 1.11  1999/09/29 09:24:20  fca
63 Introduction of the Copyright and cvs Log
64
65 */
66
67 ///////////////////////////////////////////////////////////////////////////////
68 //
69 //  Inner Traking System version 5
70 //  This class contains the base procedures for the Inner Tracking System
71 //
72 // Authors: R. Barbera, B. S. Nilsen.
73 // version 5.
74 // Created September 17 1999.
75 //
76 ///////////////////////////////////////////////////////////////////////////////
77 #include <iostream.h>
78 #include <iomanip.h>
79 #include <stdio.h>
80 #include <stdlib.h>
81 #include <TMath.h>
82 #include <TGeometry.h>
83 #include <TNode.h>
84 #include <TTUBE.h>
85 #include <TFile.h>    // only required for Tracking function?
86 #include <TCanvas.h>
87 #include <TObjArray.h>
88 #include <TLorentzVector.h>
89 #include <TObjString.h>
90 #include <TClonesArray.h>
91 #include <TBRIK.h>
92 #include <TSystem.h>
93
94 #include "AliMC.h"
95 #include "AliRun.h"
96 #include "../TGeant3/TGeant3.h"
97 #include "AliITShit.h"
98 #include "AliITSGeant3Geometry.h"
99 #include "AliITS.h"
100 #include "AliITSv5.h"
101 #include "AliITSgeom.h"
102 #include "AliITSgeomSPD.h"
103 #include "AliITSgeomSDD.h"
104 #include "AliITSgeomSSD.h"
105
106 ClassImp(AliITSv5)
107  
108 //_____________________________________________________________________________
109 AliITSv5::AliITSv5() {
110 ////////////////////////////////////////////////////////////////////////
111 //    Standard default constructor for the ITS version 5.
112 ////////////////////////////////////////////////////////////////////////
113     Int_t i;
114
115     fIdN    = 0;
116     fIdName = 0;
117     fIdSens = 0;
118     fEuclidOut    = kFALSE; // Don't write Euclide file
119     fGeomDetOut   = kFALSE; // Don't write .det file
120     fGeomDetIn    = kTRUE; // Read .det file
121     fGeomOldDetIn = kTRUE;  // Read old formatted .det file
122     fMajorVersion = IsVersion();
123     fMinorVersion = 1;
124     for(i=0;i<60;i++) fRead[i] = '\0';
125     for(i=0;i<60;i++) fWrite[i] = '\0';
126     for(i=0;i<60;i++) fEuclidGeomDet[i] = '\0';
127 }
128 //_____________________________________________________________________________
129 AliITSv5::AliITSv5(const char *name, const char *title) : AliITS(name, title){
130 ////////////////////////////////////////////////////////////////////////
131 //    Standard constructor for the ITS version 5.
132 ////////////////////////////////////////////////////////////////////////
133     Int_t i;
134
135     fIdN    = 6;
136     fIdName    = new TString[fIdN];
137     fIdName[0] = "ITS1";
138     fIdName[1] = "ITS2";
139     fIdName[2] = "ITS3";
140     fIdName[3] = "ITS4";
141     fIdName[4] = "ITS5";
142     fIdName[5] = "ITS6";
143     fIdSens    = new Int_t[fIdN];
144     for (i=0;i<fIdN;i++) fIdSens[i] = 0;
145     fEuclidOut    = kFALSE; // Don't write Euclide file
146     fGeomDetOut   = kFALSE; // Don't write .det file
147     fGeomDetIn    = kTRUE; // Read .det file
148     fGeomOldDetIn = kTRUE;  // Read old formatted .det file
149     fMajorVersion = IsVersion();
150     fMinorVersion = 1;
151     for(i=0;i<60;i++) fRead[i] = '\0';
152     for(i=0;i<60;i++) fWrite[i] = '\0';
153
154     fEuclidMaterial = "$ALICE_ROOT/Euclid/ITSgeometry_5.tme";
155     fEuclidGeometry = "$ALICE_ROOT/Euclid/ITSgeometry_5.euc";
156     strncpy(fEuclidGeomDet,"$ALICE_ROOT/ITS/ITSgeometry_v5.det",60);
157     strncpy(fRead,fEuclidGeomDet,60);
158     strncpy(fWrite,fEuclidGeomDet,60);
159 }
160 //____________________________________________________________________________
161 AliITSv5::AliITSv5(const AliITSv5 &source){
162 ////////////////////////////////////////////////////////////////////////
163 //     Copy Constructor for ITS version 5.
164 ////////////////////////////////////////////////////////////////////////
165     if(&source == this) return;
166     Warning("Copy Constructor","Not allowed to copy AliITSv5");
167     return;
168 }
169 //_____________________________________________________________________________
170 AliITSv5& AliITSv5::operator=(const AliITSv5 &source){
171 ////////////////////////////////////////////////////////////////////////
172 //    Assignment operator for the ITS version 5.
173 ////////////////////////////////////////////////////////////////////////
174
175     if(&source == this) return *this;
176     Warning("= operator","Not allowed to copy AliITSv5");
177     return *this;
178
179 }
180 //_____________________________________________________________________________
181 AliITSv5::~AliITSv5() {
182 ////////////////////////////////////////////////////////////////////////
183 //    Standard destructor for the ITS version 5.
184 ////////////////////////////////////////////////////////////////////////
185 }
186 //______________________________________________________________________
187 void AliITSv5::BuildGeometry(){
188 ////////////////////////////////////////////////////////////////////////
189 //    Geometry builder for the ITS version 5.
190 ////////////////////////////////////////////////////////////////////////
191   //
192   // Build ITS TNODE geometry for event display using detailed geometry.
193   // This function builds a simple ITS geometry used by the ROOT macro
194   // ITSdisplay.C.
195
196   TNode *top;
197   TNode *nd;
198   //const int kColorITSSPD=kRed;
199   //const int kColorITSSDD=kGreen;
200   const int kColorITSSSD=kBlue;
201   //
202   AliITSgeom  *gm = this->GetITSgeom();
203   if(gm==0) return;
204   top=gAlice->GetGeometry()->GetNode("alice");
205
206   Int_t       lay,lad,det,i;
207   Text_t      name[10];
208   Float_t     xg[3];
209   Float_t     rt[9];
210   Double_t    rtd[9];
211   TBRIK       *box;
212   TRotMatrix  *rm;
213   //TCanvas     *c1 = new TCanvas("c1","ITS");
214
215   for(lay=1;lay<=2;lay++)
216    for(lad=1;lad<=gm->GetNladders(lay);lad++)
217     for(det=1;det<=gm->GetNdetectors(lay);det++){
218           try {
219               box  = new TBRIK ("ActiveSPD","Active volume of SPD","SPD SI DET",
220                                     0.64,0.0075,4.19); 
221           } catch (...) {
222               cout << "EXCEPTION in box = new TBRIK" << endl;
223               return;
224           }
225           gm->GetTrans(lay,lad,det,xg[0],xg[1],xg[2]);
226           gm->GetRotMatrix(lay,lad,det,rt);
227           //sprintf(name,"ROT%1.1d2.2d2.2d",lay,lad,det);
228           for(i=0;i<9;i++) rtd[i] = rt[i];
229           try {
230                 rm  = new TRotMatrix(name,name,rtd);
231           } catch (...) {
232                 cout << "EXCEPTION in   new TRotMatrix" << endl;
233                 return;
234           }
235          top->cd();
236           //sprintf(name,"ND%1.1d2.2d2.2d",lay,lad,det); 
237          try {
238               nd  = new TNode("SPD"," ",box,xg[0],xg[1],xg[2],rm);
239          } catch (...) {
240               cout << "EXCEPTION in new TNode" << endl;
241               return;
242          }
243          nd->SetLineColor(kColorITSSSD);
244          fNodes->Add(nd);
245     }
246
247   for(lay=3;lay<=3;lay++)
248    for(lad=1;lad<=gm->GetNladders(lay);lad++)
249     for(det=1;det<=gm->GetNdetectors(lay);det++){
250           try {
251               box  = new TBRIK ("ActiveSDD","Active volume of SDD","SDD SI DET",
252                                     3.5,0.014,3.763); 
253           } catch (...) {
254               cout << "EXCEPTION in box = new TBRIK" << endl;
255               return;
256           }
257           gm->GetTrans(lay,lad,det,xg[0],xg[1],xg[2]);
258           gm->GetRotMatrix(lay,lad,det,rt);
259           //sprintf(name,"ROT%1.1d2.2d2.2d",lay,lad,det);
260           for(i=0;i<9;i++) rtd[i] = rt[i];
261           try {
262                 rm  = new TRotMatrix(name,name,rtd);
263           } catch (...) {
264                 cout << "EXCEPTION in   new TRotMatrix" << endl;
265                 return;
266           }
267          top->cd();
268           //sprintf(name,"ND%1.1d2.2d2.2d",lay,lad,det); 
269          try {
270               nd  = new TNode("SDD"," ",box,xg[0],xg[1],xg[2],rm);
271          } catch (...) {
272               cout << "EXCEPTION in new TNode" << endl;
273               return;
274          }
275          nd->SetLineColor(kColorITSSSD);
276          fNodes->Add(nd);
277     }
278
279   for(lay=4;lay<=4;lay++)
280    for(lad=1;lad<=gm->GetNladders(lay);lad++)
281     for(det=1;det<=gm->GetNdetectors(lay);det++){
282           try {
283               box  = new TBRIK ("ActiveSDD","Active volume of SDD","SDD SI DET",
284                                     3.5,0.014,3.763); 
285           } catch (...) {
286               cout << "EXCEPTION in box = new TBRIK" << endl;
287               return;
288           }
289           gm->GetTrans(lay,lad,det,xg[0],xg[1],xg[2]);
290           gm->GetRotMatrix(lay,lad,det,rt);
291           //sprintf(name,"ROT%1.1d2.2d2.2d",lay,lad,det);
292           for(i=0;i<9;i++) rtd[i] = rt[i];
293           try {
294                 rm  = new TRotMatrix(name,name,rtd);
295           } catch (...) {
296                 cout << "EXCEPTION in   new TRotMatrix" << endl;
297                 return;
298           }
299          top->cd();
300           //sprintf(name,"ND%1.1d2.2d2.2d",lay,lad,det); 
301          try {
302               nd  = new TNode("SDD"," ",box,xg[0],xg[1],xg[2],rm);
303          } catch (...) {
304               cout << "EXCEPTION in new TNode" << endl;
305               return;
306          }
307          nd->SetLineColor(kColorITSSSD);
308          fNodes->Add(nd);
309     }
310  for(lay=5;lay<=5;lay++)
311    for(lad=1;lad<=gm->GetNladders(lay);lad++)
312     for(det=1;det<=gm->GetNdetectors(lay);det++){
313           try {
314               box  = new TBRIK ("ActiveSSD","Active volume of SSD","SSD SI DET",
315                                     3.65,0.015,2.0); 
316           } catch (...) {
317               cout << "EXCEPTION in box = new TBRIK" << endl;
318               return;
319           }
320           gm->GetTrans(lay,lad,det,xg[0],xg[1],xg[2]);
321           gm->GetRotMatrix(lay,lad,det,rt);
322           //sprintf(name,"ROT%1.1d2.2d2.2d",lay,lad,det);
323           for(i=0;i<9;i++) rtd[i] = rt[i];
324           try {
325                 rm  = new TRotMatrix(name,name,rtd);
326           } catch (...) {
327                 cout << "EXCEPTION in   new TRotMatrix" << endl;
328                 return;
329           }
330          top->cd();
331           //sprintf(name,"ND%1.1d2.2d2.2d",lay,lad,det); 
332          try {
333               nd  = new TNode("SSD"," ",box,xg[0],xg[1],xg[2],rm);
334          } catch (...) {
335               cout << "EXCEPTION in new TNode" << endl;
336               return;
337          }
338          nd->SetLineColor(kColorITSSSD);
339          fNodes->Add(nd);
340     }
341
342  for(lay=6;lay<=6;lay++)
343    for(lad=1;lad<=gm->GetNladders(lay);lad++)
344     for(det=1;det<=gm->GetNdetectors(lay);det++){
345           try {
346               box  = new TBRIK ("ActiveSSD","Active volume of SSD","SSD SI DET",
347                                     3.65,0.015,2.0); 
348           } catch (...) {
349               cout << "EXCEPTION in box = new TBRIK" << endl;
350               return;
351           }
352
353           gm->GetTrans(lay,lad,det,xg[0],xg[1],xg[2]); 
354           gm->GetRotMatrix(lay,lad,det,rt);
355           //sprintf(name,"ROT%1.1d2.2d2.2d",lay,lad,det);
356           for(i=0;i<9;i++) rtd[i] = rt[i];
357           try {
358                 rm  = new TRotMatrix(name,name,rtd);
359           } catch (...) {
360                 cout << "EXCEPTION in   new TRotMatrix" << endl;
361                 return;
362           }
363          top->cd();
364           //sprintf(name,"ND%1.1d2.2d2.2d",lay,lad,det); 
365          try {
366               nd  = new TNode("SSD"," ",box,xg[0],xg[1],xg[2],rm);
367          } catch (...) {
368               cout << "EXCEPTION in new TNode" << endl;
369               return;
370          }
371          nd->SetLineColor(kColorITSSSD);
372          fNodes->Add(nd);
373     }
374
375
376
377 //_____________________________________________________________________________
378 void AliITSv5::CreateMaterials(){
379 ////////////////////////////////////////////////////////////////////////
380 //     Read a file containing the materials for the ITS version 5.
381 ////////////////////////////////////////////////////////////////////////
382     char *filtmp;
383
384   filtmp = gSystem->ExpandPathName(fEuclidMaterial.Data());
385
386   FILE *file = fopen(filtmp,"r");
387   if(file) {
388     fclose(file);
389     ReadEuclidMedia(filtmp);
390   } else {
391     Error("CreateMaterials"," THE MEDIA FILE %s DOES NOT EXIST !",filtmp);
392     exit(1);
393   } // end if(file)
394 }
395 //_____________________________________________________________________________
396 void AliITSv5::CreateGeometry(){
397 //////////////////////////////////////////////////////////////////////
398 //    This is the geometry used for the ITS Pre-TDR and comes from an 
399 // Euclid to Geant conversion. The only difference
400 // is in the details of the ITS supports. The detectors elements, 
401 // detector numbering, and local and global reference frames are shown in
402 //  the following figures.
403 //Begin_Html
404 /*
405 <img src="picts/ITS/its1+2_convention_front_5.gif">
406 </pre>
407 <br clear=left>
408 <font size=+2 color=red>
409 <p>This shows the front view of the SPDs.
410 </font>
411 <pre>
412 <img src="picts/ITS/its1+2_convention_side_5.gif">
413 </pre>
414 <br clear=left>
415 <font size=+2 color=red>
416 <p>This shows the perspective view of the SPDs.
417 </font>
418 <img src="picts/ITS/its1+2_tree.gif">
419 </pre>
420 <br clear=left>
421 <font size=+2 color=red>
422 <p>This shows the geometry Tree for the SPDs.
423 </font>
424 <pre>
425
426 <pre>
427 <img src="picts/ITS/its3+4_convention_front_5.gif">
428 </pre>
429 <br clear=left>
430 <font size=+2 color=red>
431 <p>This shows the front view of the SDDs.
432 </font>
433 <pre>
434 <img src="picts/ITS/its3+4_convention_side_5.gif">
435 </pre>
436 <br clear=left>
437 <font size=+2 color=red>
438 <p>This shows the perspective view of the SDDs.
439 </font>
440 <img src="picts/ITS/its3+4_tree.gif">
441 </pre>
442 <br clear=left>
443 <font size=+2 color=red>
444 <p>This shows the geometry Tree for the SDDs.
445 </font>
446 <pre>
447
448 <pre>
449 <img src="picts/ITS/its5+6_convention_front_5.gif">
450 </pre>
451 <br clear=left>
452 <font size=+2 color=red>
453 <p>This shows the front view of the SSDs.
454 </font>
455 <pre>
456 <img src="picts/ITS/its5+6_convention_side_5.gif">
457 </pre>
458 <br clear=left>
459 <font size=+2 color=red>
460 <p>This shows the perspective view of the SSDs.
461 </font>
462 <pre>
463 <img src="picts/ITS/its5+6_tree.gif">
464 </pre>
465 <br clear=left>
466 <font size=+2 color=red>
467 <p>This shows the geometry Tree for the SSDs.
468 </font>
469 <pre>
470
471
472 <img src="picts/ITS/its_layer1-6_2.gif">
473 </pre>
474 <br clear=left>
475 <font size=+2 color=red>
476 <p>This shows the front view of the whole ITS..
477 </font>
478 <pre>
479
480 <img src="picts/ITS/its_layer1-6_1.gif">
481 </pre>
482 <br clear=left>
483 <font size=+2 color=red>
484 <p>This shows the perspective view of the whole ITS..
485 </font>
486 <pre>
487
488 <img src="picts/ITS/its1-6_tree.gif">
489 </pre>
490 <br clear=left>
491 <font size=+2 color=red>
492 <p>This shows the geometry Tree for the whole ITS.
493 </font>
494 <pre>
495 */
496 //End_Html
497 //
498 //
499 //      Here are shown the details of the ITS support cones and services.
500 // First is a GEANT tree showing the organization of all of the volumes
501 // that make up the ITS supports and services.
502 //Begin_Html
503 /*
504 <img src="picts/ITS/supports_tree.gif">
505  */
506 //End_Html
507 //     What follows are a number of figures showing what these support
508 // structures look like.
509 //Begin_Html
510 /*
511
512 <img src="picts/ITS/supports_3.gif">
513 </pre>
514 <br clear=left>
515 <font size=+2 color=red>
516 <p>This shows the geometry of the supports for the Drift and Strip layers only.
517 </font>
518 <pre>
519
520 <img src="picts/ITS/supports_2.gif">
521 </pre>
522 <br clear=left>
523 <font size=+2 color=red>
524 <p>This shows the geometry of the supports for the Drift and Strip layers in front cut out.
525 </font>
526 <pre>
527
528 <img src="picts/ITS/supports_1.gif">
529 </pre>
530 <br clear=left>
531 <font size=+2 color=red>
532 <p>This shows the geometry of the supports for the Drift and Strip layers in a back cut out..
533 </font>
534 <pre>
535
536 <img src="picts/ITS/suppssd.gif">
537 </pre>
538 <br clear=left>
539 <font size=+2 color=red>
540 <p>This shows the geometry for the Strip layers supports.
541 </font>
542 <pre>
543
544 <img src="picts/ITS/suppsdd.gif">
545 </pre>
546 <br clear=left>
547 <font size=+2 color=red>
548 <p>This shows the geometry for the Drift layers supports.
549 </font>
550 <pre>
551  */
552 //End_Html
553 //
554 //    Read a file containing the geometry for the ITS version 5.
555 ////////////////////////////////////////////////////////////////////////
556
557
558     char topvol[5];
559     char *filtmp;
560
561   filtmp = gSystem->ExpandPathName(fEuclidGeometry.Data());
562   FILE *file = fopen(filtmp,"r");
563   delete [] filtmp;
564   if(file) {
565     fclose(file);
566     cout << "Ready to read Euclid geometry file" << endl;
567     ReadEuclid(fEuclidGeometry.Data(),topvol);
568     printf("Read in euclid geometries\n");
569   } else {
570     Error("CreateGeometry"," THE GEOM FILE %s DOES NOT EXIST !",
571           fEuclidGeometry.Data());
572     exit(1);
573   } // end if(file)
574   //
575   // Place the ITS ghost volume ITSV in its mother volume (ALIC) and make it
576   // invisible
577   //
578   gMC->Gspos("ITSV",1,"ALIC",0,0,0,0,"ONLY");
579   //
580   // Outputs the geometry tree in the EUCLID/CAD format if requested to do so
581   
582     if (fEuclidOut) {
583       gMC->WriteEuclid("ITSgeometry", "ITSV", 1, 5);
584     } // end if (fEuclidOut)
585
586     cout << "finished with euclid geometrys" << endl;
587 }
588 //______________________________________________________________________
589 void AliITSv5::ReadOldGeometry(const char *filename){
590     // read in the file containing the transformations for the active
591     // volumes for the ITS version 5. This is expected to be in a file
592     // ending in .det. This geometry is kept in the AliITSgeom class.
593     Int_t size;
594     char *filtmp;
595     FILE *file;
596
597     if(fITSgeom!=0) delete fITSgeom;
598     filtmp = gSystem->ExpandPathName(filename);
599     size = strlen(filtmp);
600     if(size>4 && fGeomDetIn){
601         filtmp[size-3] = 'd'; // change from .euc to .det
602         filtmp[size-2] = 'e';
603         filtmp[size-1] = 't';
604         file = fopen(filtmp,"r");
605         if(file){ // if file exists use it to fill AliITSgeom structure.
606             fclose(file);
607             fITSgeom = new AliITSgeom(filtmp);
608             fITSgeom->DefineShapes(3); // if fShape isn't defined define it.
609             // Now define the detector types/shapes.
610             fITSgeom->ReSetShape(kSPD,new AliITSgeomSPD300());
611             fITSgeom->ReSetShape(kSDD,new AliITSgeomSDD300());
612             fITSgeom->ReSetShape(kSSD,new AliITSgeomSSD175());
613         }else{
614             fITSgeom = 0;
615             // fill AliITSgeom structure from geant structure just filled above
616         }// end if(file)
617         delete [] filtmp;
618     }// end if(size>4)
619 }
620 //______________________________________________________________________
621 void AliITSv5::InitAliITSgeom(){
622 //     Based on the geometry tree defined in Geant 3.21, this
623 // routine initilizes the Class AliITSgeom from the Geant 3.21 ITS geometry
624 // sturture.
625     if(!((TGeant3*)gMC)) {
626         Error("InitAliITSgeom",
627                 "Wrong Monte Carlo. InitAliITSgeom uses TGeant3 calls");
628         return;
629     } // end if
630     cout << "Reading Geometry transformation directly from Geant 3." << endl;
631     const Int_t nlayers = 6;
632     const Int_t ndeep = 7;
633     Int_t itsGeomTreeNames[nlayers][ndeep],lnam[20],lnum[20];
634     Int_t nlad[nlayers],ndet[nlayers];
635     Double_t t[3],r[10];
636     Float_t  par[20],att[20];
637     Int_t    npar,natt,idshape,imat,imed;
638     AliITSGeant3Geometry *ig = new AliITSGeant3Geometry();
639     Int_t mod,lay,lad,det,i,j,k;
640     char *names[nlayers][ndeep] = {
641         {"ALIC","ITSV","ITSD","IT12","I132","I186","ITS1"}, // lay=1
642         {"ALIC","ITSV","ITSD","IT12","I132","I131","ITS2"}, // lay=2
643         {"ALIC","ITSV","ITSD","IT34","I004","I302","ITS3"}, // lay=3
644         {"ALIC","ITSV","ITSD","IT34","I005","I402","ITS4"}, // lay=4
645         {"ALIC","ITSV","ITSD","IT56","I565","I562","ITS5"}, // lay=5
646         {"ALIC","ITSV","ITSD","IT56","I569","I566","ITS6"}};// lay=6
647     Int_t itsGeomTreeCopys[nlayers][ndeep] = {{1,1,1,1,10, 2,4}, // lay=1
648                                               {1,1,1,1,10, 4,4}, // lay=2
649                                               {1,1,1,1,14, 6,1}, // lay=3
650                                               {1,1,1,1,22, 8,1}, // lay=4
651                                               {1,1,1,1,34,23,1}, // lay=5
652                                               {1,1,1,1,38,26,1}};// lay=6
653
654     // Sorry, but this is not very pritty code. It should be replaced
655     // at some point with a version that can search through the geometry
656     // tree its self.
657     for(i=0;i<20;i++) lnam[i] = lnum[i] = 0;
658     for(i=0;i<nlayers;i++)for(j=0;j<ndeep;j++) 
659         itsGeomTreeNames[i][j] = ig->StringToInt(names[i][j]);
660     mod = 0;
661     for(i=0;i<nlayers;i++){
662         k = 1;
663         for(j=0;j<ndeep;j++) if(itsGeomTreeCopys[i][j]!=0)
664             k *= TMath::Abs(itsGeomTreeCopys[i][j]);
665         mod += k;
666     } // end for i
667
668     if(fITSgeom!=0) delete fITSgeom;
669     nlad[0]=20;nlad[1]=40;nlad[2]=14;nlad[3]=22;nlad[4]=34;nlad[5]=38;
670     ndet[0]=4;ndet[1]=4;ndet[2]=6;ndet[3]=8;ndet[4]=22;ndet[5]=25;
671     fITSgeom = new AliITSgeom(0,6,nlad,ndet,mod);
672     mod = -1;
673     for(lay=1;lay<=nlayers;lay++){
674         for(j=0;j<ndeep;j++) lnam[j] = itsGeomTreeNames[lay-1][j];
675         for(j=0;j<ndeep;j++) lnum[j] = itsGeomTreeCopys[lay-1][j];
676         switch (lay){
677         case 1: case 2: // layers 1 and 2 are a bit special
678             lad = 0;
679             for(j=1;j<=itsGeomTreeCopys[lay-1][4];j++){
680                 lnum[4] = j;
681                 for(k=1;k<=itsGeomTreeCopys[lay-1][5];k++){
682                     lad++;
683                     lnum[5] = k;
684                     for(det=1;det<=itsGeomTreeCopys[lay-1][6];det++){
685                         lnum[6] = det;
686                         mod++;
687                         ig->GetGeometry(ndeep,lnam,lnum,t,r,idshape,npar,natt,
688                                         par,att,imat,imed);
689                         fITSgeom->CreatMatrix(mod,lay,lad,det,kSPD,t,r);
690                         if(!(fITSgeom->IsShapeDefined((Int_t)kSPD)))
691                              fITSgeom->ReSetShape(kSPD,
692                                                   new AliITSgeomSPD300());
693                     } // end for det
694                 } // end for k
695             } // end for j
696             break;
697         case 3: case 4: case 5: case 6: // layers 3-6
698             lnum[6] = 1;
699             for(lad=1;lad<=itsGeomTreeCopys[lay-1][4];lad++){
700                 lnum[4] = lad;
701                 for(det=1;det<=itsGeomTreeCopys[lay-1][5];det++){
702                     lnum[5] = det;
703                     mod++;
704                     ig->GetGeometry(7,lnam,lnum,t,r,idshape,npar,natt,
705                                     par,att,imat,imed);
706                     switch (lay){
707                     case 3: case 4:
708                         fITSgeom->CreatMatrix(mod,lay,lad,det,kSDD,t,r);
709                         if(!(fITSgeom->IsShapeDefined(kSDD))) 
710                             fITSgeom->ReSetShape(kSDD,new AliITSgeomSDD300());
711                         break;
712                     case 5: case 6:
713                         fITSgeom->CreatMatrix(mod,lay,lad,det,kSSD,t,r);
714                         if(!(fITSgeom->IsShapeDefined(kSSD))) 
715                             fITSgeom->ReSetShape(kSSD,new AliITSgeomSSD175());
716                         break;
717                         } // end switch
718                 } // end for det
719             } // end for lad
720             break;
721         } // end switch
722     } // end for lay
723     return;
724 }
725 //_____________________________________________________________________________
726 void AliITSv5::Init(){
727 ////////////////////////////////////////////////////////////////////////
728 //     Initialise the ITS after it has been created.
729 ////////////////////////////////////////////////////////////////////////
730     Int_t i;
731     Bool_t bg = kFALSE;
732
733     cout << endl;
734     for(i=0;i<30;i++) cout << "*";cout << " ITSv5_Init ";
735     for(i=0;i<30;i++) cout << "*";cout << endl;
736 //
737     if(fITSgeom==0) bg = kTRUE;
738     if(fRead[0]=='\0') strncpy(fRead,fEuclidGeomDet,60);
739     if(fWrite[0]=='\0') strncpy(fWrite,fEuclidGeomDet,60);
740     if(fGeomDetIn && !fGeomOldDetIn){
741         if(fITSgeom!=0) delete fITSgeom;
742         fITSgeom = new AliITSgeom();
743         fITSgeom->ReadNewFile(fRead);
744     } // end if
745     if(fGeomDetIn &&  fGeomOldDetIn) ReadOldGeometry(fEuclidGeometry.Data());
746
747     if(!fGeomDetIn) this->InitAliITSgeom();
748     if(fGeomDetOut) fITSgeom->WriteNewFile(fWrite);
749     if(bg) BuildGeometry(); // call BuildGeometry here if fITSgeom was not
750                             // defined ealier
751     AliITS::Init();
752 //
753     for(i=0;i<72;i++) cout << "*";
754     cout << endl;
755 }
756 //_____________________________________________________________________________
757 void AliITSv5::StepManager(){
758 ////////////////////////////////////////////////////////////////////////
759 //    Called for every step in the ITS, then calles the AliITShit class
760 // creator with the information to be recoreded about that hit.
761 //     The value of the macro ALIITSPRINTGEOM if set to 1 will allow the
762 // printing of information to a file which can be used to create a .det
763 // file read in by the routine CreateGeometry(). If set to 0 or any other
764 // value except 1, the default behavior, then no such file is created nor
765 // it the extra variables and the like used in the printing allocated.
766 ////////////////////////////////////////////////////////////////////////
767   Int_t          copy, id;
768   Int_t          copy1,copy2;
769   Float_t        hits[8];
770   Int_t          vol[4];
771   TLorentzVector position, momentum;
772   TClonesArray   &lhits = *fHits;
773   //
774   // Track status
775   vol[3] = 0;
776   if(gMC->IsTrackInside())      vol[3] +=  1;
777   if(gMC->IsTrackEntering())    vol[3] +=  2;
778   if(gMC->IsTrackExiting())     vol[3] +=  4;
779   if(gMC->IsTrackOut())         vol[3] +=  8;
780   if(gMC->IsTrackDisappeared()) vol[3] += 16;
781   if(gMC->IsTrackStop())        vol[3] += 32;
782   if(gMC->IsTrackAlive())       vol[3] += 64;
783   //
784   // Fill hit structure.
785   if(!(gMC->TrackCharge())) return;
786   //
787   // Only entering charged tracks
788   if((id = gMC->CurrentVolID(copy)) == fIdSens[0]) {
789       vol[0] = 1;
790       id = gMC->CurrentVolOffID(0,copy);
791       //detector copy in the ladder = 1<->4  (ITS1)
792       vol[1] = copy;
793       gMC->CurrentVolOffID(1,copy1);
794       //ladder copy in the module   = 1<->2  (I186)
795       gMC->CurrentVolOffID(2,copy2);
796       //module copy in the layer    = 1<->10 (I132)
797       vol[2] = copy1+(copy2-1)*2;//# of ladders in one module  = 2
798   } else if(id == fIdSens[1]){
799       vol[0] = 2;
800       id = gMC->CurrentVolOffID(0,copy);
801       //detector copy in the ladder = 1<->4  (ITS2)
802       vol[1] = copy;
803       gMC->CurrentVolOffID(1,copy1);
804       //ladder copy in the module   = 1<->4  (I131)
805       gMC->CurrentVolOffID(2,copy2);
806       //module copy in the layer    = 1<->10 (I132)
807       vol[2] = copy1+(copy2-1)*4;//# of ladders in one module  = 4
808   } else if(id == fIdSens[2]){
809       vol[0] = 3;
810       id = gMC->CurrentVolOffID(1,copy);
811       //detector copy in the ladder = 1<->5  (ITS3 is inside I314)
812       vol[1] = copy;
813       id = gMC->CurrentVolOffID(2,copy);
814       //ladder copy in the layer    = 1<->12 (I316)
815       vol[2] = copy;
816   } else if(id == fIdSens[3]){
817       vol[0] = 4;
818       id = gMC->CurrentVolOffID(1,copy);
819       //detector copy in the ladder = 1<->8  (ITS4 is inside I414)
820       vol[1] = copy;
821       id = gMC->CurrentVolOffID(2,copy);
822       //ladder copy in the layer    = 1<->22 (I417)
823       vol[2] = copy;
824   }else if(id == fIdSens[4]){
825       vol[0] = 5;
826       id = gMC->CurrentVolOffID(1,copy);
827       //detector copy in the ladder = 1<->23  (ITS5 is inside I562)
828       vol[1] = copy;
829       id = gMC->CurrentVolOffID(2,copy);
830      //ladder copy in the layer    = 1<->34 (I565)
831       vol[2] = copy;
832   }else if(id == fIdSens[5]){
833       vol[0] = 6;
834       id = gMC->CurrentVolOffID(1,copy);
835       //detector copy in the ladder = 1<->26  (ITS6 is inside I566)
836       vol[1] = copy;
837       id = gMC->CurrentVolOffID(2,copy);
838       //ladder copy in the layer = 1<->38 (I569)
839       vol[2] = copy;
840   } else {
841       return; // not an ITS volume?
842   } // end if/else if (gMC->CurentVolID(copy) == fIdSens[i])
843 //
844   gMC->TrackPosition(position);
845   gMC->TrackMomentum(momentum);
846   hits[0]=position[0];
847   hits[1]=position[1];
848   hits[2]=position[2];
849   hits[3]=momentum[0];
850   hits[4]=momentum[1];
851   hits[5]=momentum[2];
852   hits[6]=gMC->Edep();
853   hits[7]=gMC->TrackTime();
854   // Fill hit structure with this new hit.
855   new(lhits[fNhits++]) AliITShit(fIshunt,gAlice->CurrentTrack(),vol,hits);
856   return;
857 }