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