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