]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSv5asymm.cxx
Implementing ESD functionality in the NewIO (Yu.Belikov)
[u/mrichter/AliRoot.git] / ITS / AliITSv5asymm.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 with asymmetric services
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 October 7 2000.
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 "AliITSv5asymm.h"
55 #include "AliRun.h"
56
57 ClassImp(AliITSv5asymm)
58  
59 //_____________________________________________________________________________
60 AliITSv5asymm::AliITSv5asymm() {
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 = 3;
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 AliITSv5asymm::AliITSv5asymm(const char *name, const char *title) : AliITS(name, title){
81 /////////////////////////////////////////////////////////////////////////////
82 //    Standard constructor for the ITS version 5 with symmetrical services.
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 = 3;
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_5asymm.tme";
106     fEuclidGeometry = "$ALICE_ROOT/Euclid/ITSgeometry_5asymm.euc";
107     strncpy(fEuclidGeomDet,"$ALICE_ROOT/Euclid/ITSgeometry_5.det",60);
108     strncpy(fRead,fEuclidGeomDet,60);
109     strncpy(fWrite,fEuclidGeomDet,60);
110 }
111 //____________________________________________________________________________
112 AliITSv5asymm::AliITSv5asymm(const AliITSv5asymm &source){
113 ////////////////////////////////////////////////////////////////////////
114 //     Copy Constructor for ITS version 5.
115 ////////////////////////////////////////////////////////////////////////
116     if(&source == this) return;
117     Warning("Copy Constructor","Not allowed to copy AliITSv5asymm");
118     return;
119 }
120 //_____________________________________________________________________________
121 AliITSv5asymm& AliITSv5asymm::operator=(const AliITSv5asymm &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 AliITSv5asymm");
128     return *this;
129
130 }
131 //_____________________________________________________________________________
132 AliITSv5asymm::~AliITSv5asymm() {
133 ////////////////////////////////////////////////////////////////////////
134 //    Standard destructor for the ITS version 5.
135 ////////////////////////////////////////////////////////////////////////
136 }
137 //_____________________________________________________________________________
138 void AliITSv5asymm::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 AliITSv5asymm::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 AliITSv5asymm::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     char topvol[5];
508     char *filtmp;
509
510   filtmp = gSystem->ExpandPathName(fEuclidGeometry.Data());
511   FILE *file = fopen(filtmp,"r");
512   delete [] filtmp;
513   if(file) {
514     fclose(file);
515     if(fDebug) cout << ClassName() << ": Ready to read Euclid geometry file" << endl;
516     ReadEuclid(fEuclidGeometry.Data(),topvol);
517     if(fDebug) cout << ClassName() << ": Read in euclid geometries" << endl;
518   } else 
519     Fatal("CreateGeometry"," THE GEOM FILE %s DOES NOT EXIST !",
520           fEuclidGeometry.Data());
521   // end if(file)
522   //
523   // Place the ITS ghost volume ITSV in its mother volume (ALIC) and make it
524   // invisible
525   //
526   gMC->Gspos("ITSV",1,"ALIC",0,0,0,0,"ONLY");
527   //
528   // Outputs the geometry tree in the EUCLID/CAD format if requested to do so
529   
530     if (fEuclidOut) {
531       gMC->WriteEuclid("ITSgeometry", "ITSV", 1, 5);
532     } // end if (fEuclidOut)
533
534     if(fDebug) cout << ClassName() << ": finished with euclid geometrys" << endl;
535 }
536
537 //______________________________________________________________________
538 void AliITSv5asymm::ReadOldGeometry(const char *filename){
539     // read in the file containing the transformations for the active
540     // volumes for the ITS version 5. This is expected to be in a file
541     // ending in .det. This geometry is kept in the AliITSgeom class.
542     Int_t size;
543     char *filtmp;
544     FILE *file;
545
546     filtmp = gSystem->ExpandPathName(filename);
547     size = strlen(filtmp);
548     if(size>4 && fGeomDetIn){
549         filtmp[size-3] = 'd'; // change from .euc to .det
550         filtmp[size-2] = 'e';
551         filtmp[size-1] = 't';
552         file = fopen(filtmp,"r");
553         if(file){ // if file exists use it to fill AliITSgeom structure.
554             fclose(file);
555             fITSgeom = new AliITSgeom(filtmp);
556             fITSgeom->DefineShapes(3); // if fShape isn't defined define it.
557             // Now define the detector types/shapes.
558             fITSgeom->ReSetShape(kSPD,(TObject *) new AliITSgeomSPD300());
559             fITSgeom->ReSetShape(kSDD,(TObject *) new AliITSgeomSDD300());
560             fITSgeom->ReSetShape(kSSD,(TObject *) new AliITSgeomSSD());
561         }else{
562             fITSgeom = 0;
563             // fill AliITSgeom structure from geant structure just filled above
564         }// end if(file)
565         delete [] filtmp;
566     }// end if(size>4)
567 }
568 //______________________________________________________________________
569 void AliITSv5asymm::InitAliITSgeom(){
570 //     Based on the geometry tree defined in Geant 3.21, this
571 // routine initilizes the Class AliITSgeom from the Geant 3.21 ITS geometry
572 // sturture.
573 //    if(gMC->IsA()!=TGeant3::Class()) {
574     if(strcmp(gMC->GetName(),"TGeant3")) {
575         Error("InitAliITSgeom",
576                 "Wrong Monte Carlo. InitAliITSgeom uses TGeant3 calls");
577         return;
578     } // end if
579     if(fDebug) cout << ClassName() 
580                     << ": Reading Geometry transformation directly from Geant 3." << endl;
581     const Int_t nlayers = 6;
582     const Int_t ndeep = 7;
583     Int_t itsGeomTreeNames[nlayers][ndeep],lnam[20],lnum[20];
584     Int_t nlad[nlayers],ndet[nlayers];
585     Double_t t[3],r[10];
586     Float_t  par[20],att[20];
587     Int_t    npar,natt,idshape,imat,imed;
588     AliITSGeant3Geometry *ig = new AliITSGeant3Geometry();
589     Int_t mod,lay,lad,det,i,j,k;
590     char *names[nlayers][ndeep] = {
591         {"ALIC","ITSV","ITSD","IT12","I132","I186","ITS1"}, // lay=1
592         {"ALIC","ITSV","ITSD","IT12","I132","I131","ITS2"}, // lay=2
593         {"ALIC","ITSV","ITSD","IT34","I004","I302","ITS3"}, // lay=3
594         {"ALIC","ITSV","ITSD","IT34","I005","I402","ITS4"}, // lay=4
595         {"ALIC","ITSV","ITSD","IT56","I565","I562","ITS5"}, // lay=5
596         {"ALIC","ITSV","ITSD","IT56","I569","I566","ITS6"}};// lay=6
597     Int_t itsGeomTreeCopys[nlayers][ndeep] = {{1,1,1,1,10, 2,4}, // lay=1
598                                               {1,1,1,1,10, 4,4}, // lay=2
599                                               {1,1,1,1,14, 6,1}, // lay=3
600                                               {1,1,1,1,22, 8,1}, // lay=4
601                                               {1,1,1,1,34,23,1}, // lay=5
602                                               {1,1,1,1,38,26,1}};// lay=6
603
604     // Sorry, but this is not very pritty code. It should be replaced
605     // at some point with a version that can search through the geometry
606     // tree its self.
607     if(fDebug) cout << ClassName() 
608                     << ": Reading Geometry informaton from Geant3 common blocks" << endl;
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         itsGeomTreeNames[i][j] = ig->StringToInt(names[i][j]);
612     mod = 0;
613     for(i=0;i<nlayers;i++){
614         k = 1;
615         for(j=0;j<ndeep;j++) if(itsGeomTreeCopys[i][j]!=0)
616             k *= TMath::Abs(itsGeomTreeCopys[i][j]);
617         mod += k;
618     } // end for i
619
620     if(fITSgeom!=0) delete fITSgeom;
621     nlad[0]=20;nlad[1]=40;nlad[2]=14;nlad[3]=22;nlad[4]=34;nlad[5]=38;
622     ndet[0]=4;ndet[1]=4;ndet[2]=6;ndet[3]=8;ndet[4]=23;ndet[5]=26;
623     fITSgeom = new AliITSgeom(0,6,nlad,ndet,mod);
624     mod = -1;
625     for(lay=1;lay<=nlayers;lay++){
626         for(j=0;j<ndeep;j++) lnam[j] = itsGeomTreeNames[lay-1][j];
627         for(j=0;j<ndeep;j++) lnum[j] = itsGeomTreeCopys[lay-1][j];
628         switch (lay){
629         case 1: case 2: // layers 1 and 2 are a bit special
630             lad = 0;
631             for(j=1;j<=itsGeomTreeCopys[lay-1][4];j++){
632                 lnum[4] = j;
633                 for(k=1;k<=itsGeomTreeCopys[lay-1][5];k++){
634                     lad++;
635                     lnum[5] = k;
636                     for(det=1;det<=itsGeomTreeCopys[lay-1][6];det++){
637                         lnum[6] = det;
638                         mod++;
639                         ig->GetGeometry(ndeep,lnam,lnum,t,r,idshape,npar,natt,
640                                         par,att,imat,imed);
641                         fITSgeom->CreatMatrix(mod,lay,lad,det,kSPD,t,r);
642                         if(!(fITSgeom->IsShapeDefined((Int_t)kSPD)))
643                             if(fMinorVersion==1){
644                              fITSgeom->ReSetShape(kSPD,
645                                                   new AliITSgeomSPD300());
646                             } else if(fMinorVersion==2){
647                              fITSgeom->ReSetShape(kSPD,
648                                                   new AliITSgeomSPD300());
649                             }else if(fMinorVersion==3){
650                              fITSgeom->ReSetShape(kSPD,
651                                                   new AliITSgeomSPD425Long());
652                             }else{
653                              fITSgeom->ReSetShape(kSPD,
654                                                   new AliITSgeomSPD300());
655                             } // end if
656                     } // end for det
657                 } // end for k
658             } // end for j
659             break;
660         case 3: case 4: case 5: case 6: // layers 3-6
661             lnum[6] = 1;
662             for(lad=1;lad<=itsGeomTreeCopys[lay-1][4];lad++){
663                 lnum[4] = lad;
664                 for(det=1;det<=itsGeomTreeCopys[lay-1][5];det++){
665                     lnum[5] = det;
666                     mod++;
667                     ig->GetGeometry(7,lnam,lnum,t,r,idshape,npar,natt,
668                                     par,att,imat,imed);
669                     switch (lay){
670                     case 3: case 4:
671                         fITSgeom->CreatMatrix(mod,lay,lad,det,kSDD,t,r);
672                         if(!(fITSgeom->IsShapeDefined(kSDD))) 
673                             fITSgeom->ReSetShape(kSDD,new AliITSgeomSDD300());
674                             break;
675                     case 5: case 6:
676                             fITSgeom->CreatMatrix(mod,lay,lad,det,kSSD,t,r);
677                             if(!(fITSgeom->IsShapeDefined(kSSD))) 
678                                 fITSgeom->ReSetShape(kSSD,new AliITSgeomSSD175());
679                             break;
680                         } // end switch
681                 } // end for det
682             } // end for lad
683             break;
684         } // end switch
685     } // end for lay
686     return;
687 }
688 //_____________________________________________________________________________
689 void AliITSv5asymm::Init(){
690 ////////////////////////////////////////////////////////////////////////
691 //     Initialise the ITS after it has been created.
692 ////////////////////////////////////////////////////////////////////////
693     Int_t i;
694
695     if(fDebug) {
696       cout << endl << ClassName() << ": ";
697       for(i=0;i<28;i++) cout << "*";cout << " ITSv5_Init ";
698       for(i=0;i<27;i++) cout << "*";cout << endl;
699     }
700 //
701     if(fRead[0]=='\0') strncpy(fRead,fEuclidGeomDet,60);
702     if(fWrite[0]=='\0') strncpy(fWrite,fEuclidGeomDet,60);
703     if(fITSgeom!=0) delete fITSgeom;
704     fITSgeom = new AliITSgeom();
705     if(fGeomDetIn && !fGeomOldDetIn) fITSgeom->ReadNewFile(fRead);
706     if(fGeomDetIn &&  fGeomOldDetIn) this->ReadOldGeometry(fRead);
707
708     if(!fGeomDetIn) this->InitAliITSgeom();
709     if(fGeomDetOut) fITSgeom->WriteNewFile(fWrite);
710     AliITS::Init();
711 //
712     if(fDebug) {
713       cout << ClassName() << ": ";
714       for(i=0;i<72;i++) cout << "*";
715       cout << endl;
716     }
717
718 //_____________________________________________________________________________
719 void AliITSv5asymm::StepManager(){
720 ////////////////////////////////////////////////////////////////////////
721 //    Called for every step in the ITS, then calles the AliITShit class
722 // creator with the information to be recoreded about that hit.
723 //     The value of the macro ALIITSPRINTGEOM if set to 1 will allow the
724 // printing of information to a file which can be used to create a .det
725 // file read in by the routine CreateGeometry(). If set to 0 or any other
726 // value except 1, the default behavior, then no such file is created nor
727 // it the extra variables and the like used in the printing allocated.
728 ////////////////////////////////////////////////////////////////////////
729   Int_t          copy, id;
730   Int_t          copy1,copy2;
731   Float_t        hits[8];
732   Int_t          vol[4];
733   TLorentzVector position, momentum;
734   TClonesArray   &lhits = *fHits;
735   //
736   // Track status
737   vol[3] = 0;
738   if(gMC->IsTrackInside())      vol[3] +=  1;
739   if(gMC->IsTrackEntering())    vol[3] +=  2;
740   if(gMC->IsTrackExiting())     vol[3] +=  4;
741   if(gMC->IsTrackOut())         vol[3] +=  8;
742   if(gMC->IsTrackDisappeared()) vol[3] += 16;
743   if(gMC->IsTrackStop())        vol[3] += 32;
744   if(gMC->IsTrackAlive())       vol[3] += 64;
745   //
746   // Fill hit structure.
747   if(!(gMC->TrackCharge())) return;
748   //
749   // Only entering charged tracks
750   if((id = gMC->CurrentVolID(copy)) == fIdSens[0]) {
751       vol[0] = 1;
752       id = gMC->CurrentVolOffID(0,copy);
753       //detector copy in the ladder = 1<->4  (ITS1)
754       vol[1] = copy;
755       gMC->CurrentVolOffID(1,copy1);
756       //ladder copy in the module   = 1<->2  (I186)
757       gMC->CurrentVolOffID(2,copy2);
758       //module copy in the layer    = 1<->10 (I132)
759       vol[2] = copy1+(copy2-1)*2;//# of ladders in one module  = 2
760   } else if(id == fIdSens[1]){
761       vol[0] = 2;
762       id = gMC->CurrentVolOffID(0,copy);
763       //detector copy in the ladder = 1<->4  (ITS2)
764       vol[1] = copy;
765       gMC->CurrentVolOffID(1,copy1);
766       //ladder copy in the module   = 1<->4  (I131)
767       gMC->CurrentVolOffID(2,copy2);
768       //module copy in the layer    = 1<->10 (I132)
769       vol[2] = copy1+(copy2-1)*4;//# of ladders in one module  = 4
770   } else if(id == fIdSens[2]){
771       vol[0] = 3;
772       id = gMC->CurrentVolOffID(1,copy);
773       //detector copy in the ladder = 1<->5  (ITS3 is inside I314)
774       vol[1] = copy;
775       id = gMC->CurrentVolOffID(2,copy);
776       //ladder copy in the layer    = 1<->12 (I316)
777       vol[2] = copy;
778   } else if(id == fIdSens[3]){
779       vol[0] = 4;
780       id = gMC->CurrentVolOffID(1,copy);
781       //detector copy in the ladder = 1<->8  (ITS4 is inside I414)
782       vol[1] = copy;
783       id = gMC->CurrentVolOffID(2,copy);
784       //ladder copy in the layer    = 1<->22 (I417)
785       vol[2] = copy;
786   }else if(id == fIdSens[4]){
787       vol[0] = 5;
788       id = gMC->CurrentVolOffID(1,copy);
789       //detector copy in the ladder = 1<->23  (ITS5 is inside I562)
790       vol[1] = copy;
791       id = gMC->CurrentVolOffID(2,copy);
792      //ladder copy in the layer    = 1<->34 (I565)
793       vol[2] = copy;
794   }else if(id == fIdSens[5]){
795       vol[0] = 6;
796       id = gMC->CurrentVolOffID(1,copy);
797       //detector copy in the ladder = 1<->26  (ITS6 is inside I566)
798       vol[1] = copy;
799       id = gMC->CurrentVolOffID(2,copy);
800       //ladder copy in the layer = 1<->38 (I569)
801       vol[2] = copy;
802   } else {
803       return; // not an ITS volume?
804   } // end if/else if (gMC->CurentVolID(copy) == fIdSens[i])
805 //
806   gMC->TrackPosition(position);
807   gMC->TrackMomentum(momentum);
808   hits[0]=position[0];
809   hits[1]=position[1];
810   hits[2]=position[2];
811   hits[3]=momentum[0];
812   hits[4]=momentum[1];
813   hits[5]=momentum[2];
814   hits[6]=gMC->Edep();
815   hits[7]=gMC->TrackTime();
816   // Fill hit structure with this new hit.
817   new(lhits[fNhits++]) AliITShit(fIshunt,gAlice->CurrentTrack(),vol,hits);
818   return;
819 }