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