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