1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 Revision 1.32 2001/03/23 00:12:23 nilsen
19 Set Reading of AliITSgeom data from Geant3 common blocks as the default and
20 not a .det file. Removed redundent calls to BuildGeometry.
22 Revision 1.31 2001/02/13 16:53:35 nilsen
23 Fixed a but when trying to use GEANT4. Needed to replace
24 if(!((TGeant3*)gMC)) with if(!(dynamic_casst<TGeant3*>(gMC)))
25 because just casting gMC to be TGeant3* even when it realy is a TGeant3 pointer
26 did not result in a zero value. For AliITSv5asymm and AliITSv5symm, needed
27 to fix a bug in the initilizers and a bug in BuildGeometry. This is now done
28 in the same way as in AliITSv5.cxx.
30 Revision 1.30 2001/02/09 20:06:26 nilsen
31 Fixed bug in distructor. Can't distroy fixxed length arrays. Thanks Peter.
33 Revision 1.29 2001/02/09 00:05:31 nilsen
34 Added fMajor/MinorVersion variables and made other changes to better make
35 use of the new code changes in AliITSgeom related classes.
37 Revision 1.28 2001/02/02 23:57:28 nilsen
38 Added include file that are no londer included in AliITSgeom.h
40 Revision 1.27 2001/01/30 09:23:13 hristov
41 Streamers removed (R.Brun)
43 Revision 1.26 2000/11/30 11:13:11 barbera
44 Added changes suggested by Federico Carminati on nov, 30, 2000
46 Revision 1.25 2000/10/05 20:50:00 nilsen
47 Now using root generated streamers.
49 Revision 1.14.4.12 2000/10/02 16:04:03 barbera
50 Forward declarations added
52 Revision 1.22 2000/07/10 16:07:19 fca
53 Release version of ITS code
55 Revision 1.14.4.4 2000/05/19 10:10:21 nilsen
56 fix for bug with HP and Sun unix + fix for event display in ITS-working branch
58 Revision 1.14.4.3 2000/03/04 23:46:38 nilsen
59 Fixed up the comments/documentation.
61 Revision 1.14.4.2 2000/03/02 21:53:02 nilsen
62 To make it compatable with the changes in AliRun/AliModule.
64 Revision 1.14.4.1 2000/01/12 19:03:33 nilsen
65 This is the version of the files after the merging done in December 1999.
66 See the ReadMe110100.txt file for details
68 Revision 1.14 1999/10/22 08:16:49 fca
69 Correct destructors, thanks to I.Hrivnacova
71 Revision 1.13 1999/10/06 10:15:19 fca
72 Correct bug in allocation of layer name and add destructor
74 Revision 1.12 1999/10/05 08:05:09 fca
75 Minor corrections for uninitialised variables.
77 Revision 1.11 1999/09/29 09:24:20 fca
78 Introduction of the Copyright and cvs Log
82 ///////////////////////////////////////////////////////////////////////////////
84 // Inner Traking System version 5
85 // This class contains the base procedures for the Inner Tracking System
87 // Authors: R. Barbera, B. S. Nilsen.
89 // Created September 17 1999.
91 ///////////////////////////////////////////////////////////////////////////////
97 #include <TGeometry.h>
100 #include <TFile.h> // only required for Tracking function?
102 #include <TObjArray.h>
103 #include <TLorentzVector.h>
104 #include <TObjString.h>
105 #include <TClonesArray.h>
111 #include "../TGeant3/TGeant3.h"
112 #include "AliITShit.h"
113 #include "AliITSGeant3Geometry.h"
115 #include "AliITSv5.h"
116 #include "AliITSgeom.h"
117 #include "AliITSgeomSPD.h"
118 #include "AliITSgeomSDD.h"
119 #include "AliITSgeomSSD.h"
123 //_____________________________________________________________________________
124 AliITSv5::AliITSv5() {
125 ////////////////////////////////////////////////////////////////////////
126 // Standard default constructor for the ITS version 5.
127 ////////////////////////////////////////////////////////////////////////
133 fEuclidOut = kFALSE; // Don't write Euclide file
134 fGeomDetOut = kFALSE; // Don't write .det file
135 fGeomDetIn = kFALSE; // Don't Read .det file
136 fGeomOldDetIn = kFALSE; // Don't Read old formatted .det file
137 fMajorVersion = IsVersion();
139 for(i=0;i<60;i++) fRead[i] = '\0';
140 for(i=0;i<60;i++) fWrite[i] = '\0';
141 for(i=0;i<60;i++) fEuclidGeomDet[i] = '\0';
143 //_____________________________________________________________________________
144 AliITSv5::AliITSv5(const char *name, const char *title) : AliITS(name, title){
145 ////////////////////////////////////////////////////////////////////////
146 // Standard constructor for the ITS version 5.
147 ////////////////////////////////////////////////////////////////////////
151 fIdName = new TString[fIdN];
158 fIdSens = new Int_t[fIdN];
159 for (i=0;i<fIdN;i++) fIdSens[i] = 0;
160 fEuclidOut = kFALSE; // Don't write Euclide file
161 fGeomDetOut = kFALSE; // Don't write .det file
162 fGeomDetIn = kFALSE; // Don't Read .det file
163 fGeomOldDetIn = kFALSE; // Don't Read old formatted .det file
164 fMajorVersion = IsVersion();
166 for(i=0;i<60;i++) fRead[i] = '\0';
167 for(i=0;i<60;i++) fWrite[i] = '\0';
169 fEuclidMaterial = "$ALICE_ROOT/Euclid/ITSgeometry_5.tme";
170 fEuclidGeometry = "$ALICE_ROOT/Euclid/ITSgeometry_5.euc";
171 strncpy(fEuclidGeomDet,"$ALICE_ROOT/ITS/ITSgeometry_v5.det",60);
172 strncpy(fRead,fEuclidGeomDet,60);
173 strncpy(fWrite,fEuclidGeomDet,60);
175 //____________________________________________________________________________
176 AliITSv5::AliITSv5(const AliITSv5 &source){
177 ////////////////////////////////////////////////////////////////////////
178 // Copy Constructor for ITS version 5.
179 ////////////////////////////////////////////////////////////////////////
180 if(&source == this) return;
181 Warning("Copy Constructor","Not allowed to copy AliITSv5");
184 //_____________________________________________________________________________
185 AliITSv5& AliITSv5::operator=(const AliITSv5 &source){
186 ////////////////////////////////////////////////////////////////////////
187 // Assignment operator for the ITS version 5.
188 ////////////////////////////////////////////////////////////////////////
190 if(&source == this) return *this;
191 Warning("= operator","Not allowed to copy AliITSv5");
195 //_____________________________________________________________________________
196 AliITSv5::~AliITSv5() {
197 ////////////////////////////////////////////////////////////////////////
198 // Standard destructor for the ITS version 5.
199 ////////////////////////////////////////////////////////////////////////
201 //______________________________________________________________________
202 void AliITSv5::BuildGeometry(){
203 ////////////////////////////////////////////////////////////////////////
204 // Geometry builder for the ITS version 5.
205 ////////////////////////////////////////////////////////////////////////
207 // Build ITS TNODE geometry for event display using detailed geometry.
208 // This function builds a simple ITS geometry used by the ROOT macro
213 //const int kColorITSSPD=kRed;
214 //const int kColorITSSDD=kGreen;
215 const int kColorITSSSD=kBlue;
217 AliITSgeom *gm = this->GetITSgeom();
219 top=gAlice->GetGeometry()->GetNode("alice");
228 //TCanvas *c1 = new TCanvas("c1","ITS");
230 for(lay=1;lay<=2;lay++)
231 for(lad=1;lad<=gm->GetNladders(lay);lad++)
232 for(det=1;det<=gm->GetNdetectors(lay);det++){
234 box = new TBRIK ("ActiveSPD","Active volume of SPD","SPD SI DET",
237 cout << "EXCEPTION in box = new TBRIK" << endl;
240 gm->GetTrans(lay,lad,det,xg[0],xg[1],xg[2]);
241 gm->GetRotMatrix(lay,lad,det,rt);
242 //sprintf(name,"ROT%1.1d2.2d2.2d",lay,lad,det);
243 for(i=0;i<9;i++) rtd[i] = rt[i];
245 rm = new TRotMatrix(name,name,rtd);
247 cout << "EXCEPTION in new TRotMatrix" << endl;
251 //sprintf(name,"ND%1.1d2.2d2.2d",lay,lad,det);
253 nd = new TNode("SPD"," ",box,xg[0],xg[1],xg[2],rm);
255 cout << "EXCEPTION in new TNode" << endl;
258 nd->SetLineColor(kColorITSSSD);
262 for(lay=3;lay<=3;lay++)
263 for(lad=1;lad<=gm->GetNladders(lay);lad++)
264 for(det=1;det<=gm->GetNdetectors(lay);det++){
266 box = new TBRIK ("ActiveSDD","Active volume of SDD","SDD SI DET",
269 cout << "EXCEPTION in box = new TBRIK" << endl;
272 gm->GetTrans(lay,lad,det,xg[0],xg[1],xg[2]);
273 gm->GetRotMatrix(lay,lad,det,rt);
274 //sprintf(name,"ROT%1.1d2.2d2.2d",lay,lad,det);
275 for(i=0;i<9;i++) rtd[i] = rt[i];
277 rm = new TRotMatrix(name,name,rtd);
279 cout << "EXCEPTION in new TRotMatrix" << endl;
283 //sprintf(name,"ND%1.1d2.2d2.2d",lay,lad,det);
285 nd = new TNode("SDD"," ",box,xg[0],xg[1],xg[2],rm);
287 cout << "EXCEPTION in new TNode" << endl;
290 nd->SetLineColor(kColorITSSSD);
294 for(lay=4;lay<=4;lay++)
295 for(lad=1;lad<=gm->GetNladders(lay);lad++)
296 for(det=1;det<=gm->GetNdetectors(lay);det++){
298 box = new TBRIK ("ActiveSDD","Active volume of SDD","SDD SI DET",
301 cout << "EXCEPTION in box = new TBRIK" << endl;
304 gm->GetTrans(lay,lad,det,xg[0],xg[1],xg[2]);
305 gm->GetRotMatrix(lay,lad,det,rt);
306 //sprintf(name,"ROT%1.1d2.2d2.2d",lay,lad,det);
307 for(i=0;i<9;i++) rtd[i] = rt[i];
309 rm = new TRotMatrix(name,name,rtd);
311 cout << "EXCEPTION in new TRotMatrix" << endl;
315 //sprintf(name,"ND%1.1d2.2d2.2d",lay,lad,det);
317 nd = new TNode("SDD"," ",box,xg[0],xg[1],xg[2],rm);
319 cout << "EXCEPTION in new TNode" << endl;
322 nd->SetLineColor(kColorITSSSD);
325 for(lay=5;lay<=5;lay++)
326 for(lad=1;lad<=gm->GetNladders(lay);lad++)
327 for(det=1;det<=gm->GetNdetectors(lay);det++){
329 box = new TBRIK ("ActiveSSD","Active volume of SSD","SSD SI DET",
332 cout << "EXCEPTION in box = new TBRIK" << endl;
335 gm->GetTrans(lay,lad,det,xg[0],xg[1],xg[2]);
336 gm->GetRotMatrix(lay,lad,det,rt);
337 //sprintf(name,"ROT%1.1d2.2d2.2d",lay,lad,det);
338 for(i=0;i<9;i++) rtd[i] = rt[i];
340 rm = new TRotMatrix(name,name,rtd);
342 cout << "EXCEPTION in new TRotMatrix" << endl;
346 //sprintf(name,"ND%1.1d2.2d2.2d",lay,lad,det);
348 nd = new TNode("SSD"," ",box,xg[0],xg[1],xg[2],rm);
350 cout << "EXCEPTION in new TNode" << endl;
353 nd->SetLineColor(kColorITSSSD);
357 for(lay=6;lay<=6;lay++)
358 for(lad=1;lad<=gm->GetNladders(lay);lad++)
359 for(det=1;det<=gm->GetNdetectors(lay);det++){
361 box = new TBRIK ("ActiveSSD","Active volume of SSD","SSD SI DET",
364 cout << "EXCEPTION in box = new TBRIK" << endl;
368 gm->GetTrans(lay,lad,det,xg[0],xg[1],xg[2]);
369 gm->GetRotMatrix(lay,lad,det,rt);
370 //sprintf(name,"ROT%1.1d2.2d2.2d",lay,lad,det);
371 for(i=0;i<9;i++) rtd[i] = rt[i];
373 rm = new TRotMatrix(name,name,rtd);
375 cout << "EXCEPTION in new TRotMatrix" << endl;
379 //sprintf(name,"ND%1.1d2.2d2.2d",lay,lad,det);
381 nd = new TNode("SSD"," ",box,xg[0],xg[1],xg[2],rm);
383 cout << "EXCEPTION in new TNode" << endl;
386 nd->SetLineColor(kColorITSSSD);
392 //_____________________________________________________________________________
393 void AliITSv5::CreateMaterials(){
394 ////////////////////////////////////////////////////////////////////////
395 // Read a file containing the materials for the ITS version 5.
396 ////////////////////////////////////////////////////////////////////////
399 filtmp = gSystem->ExpandPathName(fEuclidMaterial.Data());
401 FILE *file = fopen(filtmp,"r");
404 ReadEuclidMedia(filtmp);
406 Error("CreateMaterials"," THE MEDIA FILE %s DOES NOT EXIST !",filtmp);
410 //_____________________________________________________________________________
411 void AliITSv5::CreateGeometry(){
412 //////////////////////////////////////////////////////////////////////
413 // This is the geometry used for the ITS Pre-TDR and comes from an
414 // Euclid to Geant conversion. The only difference
415 // is in the details of the ITS supports. The detectors elements,
416 // detector numbering, and local and global reference frames are shown in
417 // the following figures.
420 <img src="picts/ITS/its1+2_convention_front_5.gif">
423 <font size=+2 color=red>
424 <p>This shows the front view of the SPDs.
427 <img src="picts/ITS/its1+2_convention_side_5.gif">
430 <font size=+2 color=red>
431 <p>This shows the perspective view of the SPDs.
433 <img src="picts/ITS/its1+2_tree.gif">
436 <font size=+2 color=red>
437 <p>This shows the geometry Tree for the SPDs.
442 <img src="picts/ITS/its3+4_convention_front_5.gif">
445 <font size=+2 color=red>
446 <p>This shows the front view of the SDDs.
449 <img src="picts/ITS/its3+4_convention_side_5.gif">
452 <font size=+2 color=red>
453 <p>This shows the perspective view of the SDDs.
455 <img src="picts/ITS/its3+4_tree.gif">
458 <font size=+2 color=red>
459 <p>This shows the geometry Tree for the SDDs.
464 <img src="picts/ITS/its5+6_convention_front_5.gif">
467 <font size=+2 color=red>
468 <p>This shows the front view of the SSDs.
471 <img src="picts/ITS/its5+6_convention_side_5.gif">
474 <font size=+2 color=red>
475 <p>This shows the perspective view of the SSDs.
478 <img src="picts/ITS/its5+6_tree.gif">
481 <font size=+2 color=red>
482 <p>This shows the geometry Tree for the SSDs.
487 <img src="picts/ITS/its_layer1-6_2.gif">
490 <font size=+2 color=red>
491 <p>This shows the front view of the whole ITS..
495 <img src="picts/ITS/its_layer1-6_1.gif">
498 <font size=+2 color=red>
499 <p>This shows the perspective view of the whole ITS..
503 <img src="picts/ITS/its1-6_tree.gif">
506 <font size=+2 color=red>
507 <p>This shows the geometry Tree for the whole ITS.
514 // Here are shown the details of the ITS support cones and services.
515 // First is a GEANT tree showing the organization of all of the volumes
516 // that make up the ITS supports and services.
519 <img src="picts/ITS/supports_tree.gif">
522 // What follows are a number of figures showing what these support
523 // structures look like.
527 <img src="picts/ITS/supports_3.gif">
530 <font size=+2 color=red>
531 <p>This shows the geometry of the supports for the Drift and Strip layers only.
535 <img src="picts/ITS/supports_2.gif">
538 <font size=+2 color=red>
539 <p>This shows the geometry of the supports for the Drift and Strip layers in front cut out.
543 <img src="picts/ITS/supports_1.gif">
546 <font size=+2 color=red>
547 <p>This shows the geometry of the supports for the Drift and Strip layers in a back cut out..
551 <img src="picts/ITS/suppssd.gif">
554 <font size=+2 color=red>
555 <p>This shows the geometry for the Strip layers supports.
559 <img src="picts/ITS/suppsdd.gif">
562 <font size=+2 color=red>
563 <p>This shows the geometry for the Drift layers supports.
569 // Read a file containing the geometry for the ITS version 5.
570 ////////////////////////////////////////////////////////////////////////
576 filtmp = gSystem->ExpandPathName(fEuclidGeometry.Data());
577 FILE *file = fopen(filtmp,"r");
581 cout << "Ready to read Euclid geometry file" << endl;
582 ReadEuclid(fEuclidGeometry.Data(),topvol);
583 printf("Read in euclid geometries\n");
585 Error("CreateGeometry"," THE GEOM FILE %s DOES NOT EXIST !",
586 fEuclidGeometry.Data());
590 // Place the ITS ghost volume ITSV in its mother volume (ALIC) and make it
593 gMC->Gspos("ITSV",1,"ALIC",0,0,0,0,"ONLY");
595 // Outputs the geometry tree in the EUCLID/CAD format if requested to do so
598 gMC->WriteEuclid("ITSgeometry", "ITSV", 1, 5);
599 } // end if (fEuclidOut)
601 cout << "finished with euclid geometrys" << endl;
603 //______________________________________________________________________
604 void AliITSv5::ReadOldGeometry(const char *filename){
605 // read in the file containing the transformations for the active
606 // volumes for the ITS version 5. This is expected to be in a file
607 // ending in .det. This geometry is kept in the AliITSgeom class.
612 if(fITSgeom!=0) delete fITSgeom;
613 filtmp = gSystem->ExpandPathName(filename);
614 size = strlen(filtmp);
615 if(size>4 && fGeomDetIn){
616 filtmp[size-3] = 'd'; // change from .euc to .det
617 filtmp[size-2] = 'e';
618 filtmp[size-1] = 't';
619 file = fopen(filtmp,"r");
620 if(file){ // if file exists use it to fill AliITSgeom structure.
622 fITSgeom = new AliITSgeom(filtmp);
623 fITSgeom->DefineShapes(3); // if fShape isn't defined define it.
624 // Now define the detector types/shapes.
625 fITSgeom->ReSetShape(kSPD,new AliITSgeomSPD300());
626 fITSgeom->ReSetShape(kSDD,new AliITSgeomSDD300());
627 fITSgeom->ReSetShape(kSSD,new AliITSgeomSSD175());
630 // fill AliITSgeom structure from geant structure just filled above
635 //______________________________________________________________________
636 void AliITSv5::InitAliITSgeom(){
637 // Based on the geometry tree defined in Geant 3.21, this
638 // routine initilizes the Class AliITSgeom from the Geant 3.21 ITS geometry
640 if(gMC->IsA()!=TGeant3::Class()) {
641 Error("InitAliITSgeom",
642 "Wrong Monte Carlo. InitAliITSgeom uses TGeant3 calls");
645 cout << "Reading Geometry transformation directly from Geant 3." << endl;
646 const Int_t nlayers = 6;
647 const Int_t ndeep = 7;
648 Int_t itsGeomTreeNames[nlayers][ndeep],lnam[20],lnum[20];
649 Int_t nlad[nlayers],ndet[nlayers];
651 Float_t par[20],att[20];
652 Int_t npar,natt,idshape,imat,imed;
653 AliITSGeant3Geometry *ig = new AliITSGeant3Geometry();
654 Int_t mod,lay,lad,det,i,j,k;
655 char *names[nlayers][ndeep] = {
656 {"ALIC","ITSV","ITSD","IT12","I132","I186","ITS1"}, // lay=1
657 {"ALIC","ITSV","ITSD","IT12","I132","I131","ITS2"}, // lay=2
658 {"ALIC","ITSV","ITSD","IT34","I004","I302","ITS3"}, // lay=3
659 {"ALIC","ITSV","ITSD","IT34","I005","I402","ITS4"}, // lay=4
660 {"ALIC","ITSV","ITSD","IT56","I565","I562","ITS5"}, // lay=5
661 {"ALIC","ITSV","ITSD","IT56","I569","I566","ITS6"}};// lay=6
662 Int_t itsGeomTreeCopys[nlayers][ndeep] = {{1,1,1,1,10, 2,4}, // lay=1
663 {1,1,1,1,10, 4,4}, // lay=2
664 {1,1,1,1,14, 6,1}, // lay=3
665 {1,1,1,1,22, 8,1}, // lay=4
666 {1,1,1,1,34,23,1}, // lay=5
667 {1,1,1,1,38,26,1}};// lay=6
669 // Sorry, but this is not very pritty code. It should be replaced
670 // at some point with a version that can search through the geometry
672 for(i=0;i<20;i++) lnam[i] = lnum[i] = 0;
673 for(i=0;i<nlayers;i++)for(j=0;j<ndeep;j++)
674 itsGeomTreeNames[i][j] = ig->StringToInt(names[i][j]);
676 for(i=0;i<nlayers;i++){
678 for(j=0;j<ndeep;j++) if(itsGeomTreeCopys[i][j]!=0)
679 k *= TMath::Abs(itsGeomTreeCopys[i][j]);
683 if(fITSgeom!=0) delete fITSgeom;
684 nlad[0]=20;nlad[1]=40;nlad[2]=14;nlad[3]=22;nlad[4]=34;nlad[5]=38;
685 ndet[0]=4;ndet[1]=4;ndet[2]=6;ndet[3]=8;ndet[4]=22;ndet[5]=25;
686 fITSgeom = new AliITSgeom(0,6,nlad,ndet,mod);
688 for(lay=1;lay<=nlayers;lay++){
689 for(j=0;j<ndeep;j++) lnam[j] = itsGeomTreeNames[lay-1][j];
690 for(j=0;j<ndeep;j++) lnum[j] = itsGeomTreeCopys[lay-1][j];
692 case 1: case 2: // layers 1 and 2 are a bit special
694 for(j=1;j<=itsGeomTreeCopys[lay-1][4];j++){
696 for(k=1;k<=itsGeomTreeCopys[lay-1][5];k++){
699 for(det=1;det<=itsGeomTreeCopys[lay-1][6];det++){
702 ig->GetGeometry(ndeep,lnam,lnum,t,r,idshape,npar,natt,
704 fITSgeom->CreatMatrix(mod,lay,lad,det,kSPD,t,r);
705 if(!(fITSgeom->IsShapeDefined((Int_t)kSPD)))
706 fITSgeom->ReSetShape(kSPD,
707 new AliITSgeomSPD300());
712 case 3: case 4: case 5: case 6: // layers 3-6
714 for(lad=1;lad<=itsGeomTreeCopys[lay-1][4];lad++){
716 for(det=1;det<=itsGeomTreeCopys[lay-1][5];det++){
719 ig->GetGeometry(7,lnam,lnum,t,r,idshape,npar,natt,
723 fITSgeom->CreatMatrix(mod,lay,lad,det,kSDD,t,r);
724 if(!(fITSgeom->IsShapeDefined(kSDD)))
725 fITSgeom->ReSetShape(kSDD,new AliITSgeomSDD300());
728 fITSgeom->CreatMatrix(mod,lay,lad,det,kSSD,t,r);
729 if(!(fITSgeom->IsShapeDefined(kSSD)))
730 fITSgeom->ReSetShape(kSSD,new AliITSgeomSSD175());
740 //_____________________________________________________________________________
741 void AliITSv5::Init(){
742 ////////////////////////////////////////////////////////////////////////
743 // Initialise the ITS after it has been created.
744 ////////////////////////////////////////////////////////////////////////
748 for(i=0;i<30;i++) cout << "*";cout << " ITSv5_Init ";
749 for(i=0;i<30;i++) cout << "*";cout << endl;
751 if(fRead[0]=='\0') strncpy(fRead,fEuclidGeomDet,60);
752 if(fWrite[0]=='\0') strncpy(fWrite,fEuclidGeomDet,60);
753 if(fGeomDetIn && !fGeomOldDetIn){
754 if(fITSgeom!=0) delete fITSgeom;
755 fITSgeom = new AliITSgeom();
756 fITSgeom->ReadNewFile(fRead);
758 if(fGeomDetIn && fGeomOldDetIn) ReadOldGeometry(fEuclidGeometry.Data());
760 if(!fGeomDetIn) this->InitAliITSgeom();
761 if(fGeomDetOut) fITSgeom->WriteNewFile(fWrite);
764 for(i=0;i<72;i++) cout << "*";
767 //_____________________________________________________________________________
768 void AliITSv5::StepManager(){
769 ////////////////////////////////////////////////////////////////////////
770 // Called for every step in the ITS, then calles the AliITShit class
771 // creator with the information to be recoreded about that hit.
772 // The value of the macro ALIITSPRINTGEOM if set to 1 will allow the
773 // printing of information to a file which can be used to create a .det
774 // file read in by the routine CreateGeometry(). If set to 0 or any other
775 // value except 1, the default behavior, then no such file is created nor
776 // it the extra variables and the like used in the printing allocated.
777 ////////////////////////////////////////////////////////////////////////
782 TLorentzVector position, momentum;
783 TClonesArray &lhits = *fHits;
787 if(gMC->IsTrackInside()) vol[3] += 1;
788 if(gMC->IsTrackEntering()) vol[3] += 2;
789 if(gMC->IsTrackExiting()) vol[3] += 4;
790 if(gMC->IsTrackOut()) vol[3] += 8;
791 if(gMC->IsTrackDisappeared()) vol[3] += 16;
792 if(gMC->IsTrackStop()) vol[3] += 32;
793 if(gMC->IsTrackAlive()) vol[3] += 64;
795 // Fill hit structure.
796 if(!(gMC->TrackCharge())) return;
798 // Only entering charged tracks
799 if((id = gMC->CurrentVolID(copy)) == fIdSens[0]) {
801 id = gMC->CurrentVolOffID(0,copy);
802 //detector copy in the ladder = 1<->4 (ITS1)
804 gMC->CurrentVolOffID(1,copy1);
805 //ladder copy in the module = 1<->2 (I186)
806 gMC->CurrentVolOffID(2,copy2);
807 //module copy in the layer = 1<->10 (I132)
808 vol[2] = copy1+(copy2-1)*2;//# of ladders in one module = 2
809 } else if(id == fIdSens[1]){
811 id = gMC->CurrentVolOffID(0,copy);
812 //detector copy in the ladder = 1<->4 (ITS2)
814 gMC->CurrentVolOffID(1,copy1);
815 //ladder copy in the module = 1<->4 (I131)
816 gMC->CurrentVolOffID(2,copy2);
817 //module copy in the layer = 1<->10 (I132)
818 vol[2] = copy1+(copy2-1)*4;//# of ladders in one module = 4
819 } else if(id == fIdSens[2]){
821 id = gMC->CurrentVolOffID(1,copy);
822 //detector copy in the ladder = 1<->5 (ITS3 is inside I314)
824 id = gMC->CurrentVolOffID(2,copy);
825 //ladder copy in the layer = 1<->12 (I316)
827 } else if(id == fIdSens[3]){
829 id = gMC->CurrentVolOffID(1,copy);
830 //detector copy in the ladder = 1<->8 (ITS4 is inside I414)
832 id = gMC->CurrentVolOffID(2,copy);
833 //ladder copy in the layer = 1<->22 (I417)
835 }else if(id == fIdSens[4]){
837 id = gMC->CurrentVolOffID(1,copy);
838 //detector copy in the ladder = 1<->23 (ITS5 is inside I562)
840 id = gMC->CurrentVolOffID(2,copy);
841 //ladder copy in the layer = 1<->34 (I565)
843 }else if(id == fIdSens[5]){
845 id = gMC->CurrentVolOffID(1,copy);
846 //detector copy in the ladder = 1<->26 (ITS6 is inside I566)
848 id = gMC->CurrentVolOffID(2,copy);
849 //ladder copy in the layer = 1<->38 (I569)
852 return; // not an ITS volume?
853 } // end if/else if (gMC->CurentVolID(copy) == fIdSens[i])
855 gMC->TrackPosition(position);
856 gMC->TrackMomentum(momentum);
864 hits[7]=gMC->TrackTime();
865 // Fill hit structure with this new hit.
866 new(lhits[fNhits++]) AliITShit(fIshunt,gAlice->CurrentTrack(),vol,hits);