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