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