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