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