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.16 2000/03/05 00:11:03 nilsen
21 Revision 1.15 2000/02/23 16:25:21 fca
22 AliVMC and AliGeant3 classes introduced
23 ReadEuclid moved from AliRun to AliModule
25 Revision 1.14.4.3 2000/03/04 23:46:38 nilsen
26 Fixed up the comments/documentation.
28 Revision 1.14.4.2 2000/03/02 21:53:02 nilsen
29 To make it compatable with the changes in AliRun/AliModule.
31 Revision 1.14.4.1 2000/01/12 19:03:33 nilsen
32 This is the version of the files after the merging done in December 1999.
33 See the ReadMe110100.txt file for details
35 Revision 1.14 1999/10/22 08:16:49 fca
36 Correct destructors, thanks to I.Hrivnacova
38 Revision 1.13 1999/10/06 10:15:19 fca
39 Correct bug in allocation of layer name and add destructor
41 Revision 1.12 1999/10/05 08:05:09 fca
42 Minor corrections for uninitialised variables.
44 Revision 1.11 1999/09/29 09:24:20 fca
45 Introduction of the Copyright and cvs Log
49 ///////////////////////////////////////////////////////////////////////////////
51 // Inner Traking System version 5
52 // This class contains the base procedures for the Inner Tracking System
54 // Authors: R. Barbera, B. S. Nilsen.
56 // Created September 17 1999.
58 ///////////////////////////////////////////////////////////////////////////////
60 // See AliITSv5::StepManager().
61 #define ALIITSPRINTGEOM 0 // default. don't print out gemetry information
62 //#define ALIITSPRINTGEOM 1 // print out geometry information
67 #include <TGeometry.h>
73 #if ALIITSPRINTGEOM==1
74 #include "../TGeant3/TGeant3.h"
76 #include "AliITShit.h"
79 #include "AliITSgeom.h"
83 //_____________________________________________________________________________
84 AliITSv5::AliITSv5() {
85 ////////////////////////////////////////////////////////////////////////
86 // Standard default constructor for the ITS version 5.
87 ////////////////////////////////////////////////////////////////////////
89 fId5Name = new char*[fId5N];
97 //_____________________________________________________________________________
98 AliITSv5::~AliITSv5() {
99 ////////////////////////////////////////////////////////////////////////
100 // Standard destructor for the ITS version 5.
101 ////////////////////////////////////////////////////////////////////////
104 //_____________________________________________________________________________
105 AliITSv5::AliITSv5(const char *name, const char *title) : AliITS(name, title){
106 ////////////////////////////////////////////////////////////////////////
107 // Standard constructor for the ITS version 5.
108 ////////////////////////////////////////////////////////////////////////
110 fId5Name = new char*[fId5N];
111 fId5Name[0] = "ITS1";
112 fId5Name[1] = "ITS2";
113 fId5Name[2] = "ITS3";
114 fId5Name[3] = "ITS4";
115 fId5Name[4] = "ITS5";
116 fId5Name[5] = "ITS6";
118 fEuclidMaterial = "$ALICE_ROOT/Euclid/ITSgeometry_5.tme";
119 fEuclidGeometry = "$ALICE_ROOT/Euclid/ITSgeometry_5.euc";
121 //_____________________________________________________________________________
122 void AliITSv5::BuildGeometry(){
124 // Build ITS TNODE geometry for event display using detailed geometry.
125 // This function builds a simple ITS geometry used by the ROOT macro
130 //const int kColorITS_SPD=kRed;
131 //const int kColorITS_SDD=kGreen;
132 const int kColorITS_SSD=kBlue;
134 Top=gAlice->GetGeometry()->GetNode("alice");
135 AliITSgeom *gm = this->GetITSgeom();
143 //TCanvas *c1 = new TCanvas("c1","ITS");
145 for(Int_t lay=1;lay<=2;lay++)
146 for(Int_t lad=1;lad<=gm->GetNladders(lay);lad++)
147 for(Int_t det=1;det<=gm->GetNdetectors(lay);det++){
149 box = new TBRIK ("ActiveSPD","Active volume of SPD","SPD SI DET",
152 cout << "EXCEPTION in box = new TBRIK" << endl;
155 gm->GetTrans(lay,lad,det,xg[0],xg[1],xg[2]);
156 gm->GetRotMatrix(lay,lad,det,rt);
157 //sprintf(name,"ROT%1.1d2.2d2.2d",lay,lad,det);
158 for(Int_t i=0;i<9;i++) rtd[i] = rt[i];
160 rm = new TRotMatrix(name,name,rtd);
162 cout << "EXCEPTION in new TRotMatrix" << endl;
166 //sprintf(name,"ND%1.1d2.2d2.2d",lay,lad,det);
168 nd = new TNode("SPD"," ",box,xg[0],xg[1],xg[2],rm);
170 cout << "EXCEPTION in new TNode" << endl;
173 nd->SetLineColor(kColorITS_SSD);
177 for(Int_t lay=3;lay<=3;lay++)
178 for(Int_t lad=1;lad<=gm->GetNladders(lay);lad++)
179 for(Int_t det=1;det<=gm->GetNdetectors(lay);det++){
181 box = new TBRIK ("ActiveSDD","Active volume of SDD","SDD SI DET",
184 cout << "EXCEPTION in box = new TBRIK" << endl;
187 gm->GetTrans(lay,lad,det,xg[0],xg[1],xg[2]);
188 gm->GetRotMatrix(lay,lad,det,rt);
189 //sprintf(name,"ROT%1.1d2.2d2.2d",lay,lad,det);
190 for(Int_t i=0;i<9;i++) rtd[i] = rt[i];
192 rm = new TRotMatrix(name,name,rtd);
194 cout << "EXCEPTION in new TRotMatrix" << endl;
198 //sprintf(name,"ND%1.1d2.2d2.2d",lay,lad,det);
200 nd = new TNode("SDD"," ",box,xg[0],xg[1],xg[2],rm);
202 cout << "EXCEPTION in new TNode" << endl;
205 nd->SetLineColor(kColorITS_SSD);
209 for(Int_t lay=4;lay<=4;lay++)
210 for(Int_t lad=1;lad<=gm->GetNladders(lay);lad++)
211 for(Int_t det=1;det<=gm->GetNdetectors(lay);det++){
213 box = new TBRIK ("ActiveSDD","Active volume of SDD","SDD SI DET",
216 cout << "EXCEPTION in box = new TBRIK" << endl;
219 gm->GetTrans(lay,lad,det,xg[0],xg[1],xg[2]);
220 gm->GetRotMatrix(lay,lad,det,rt);
221 //sprintf(name,"ROT%1.1d2.2d2.2d",lay,lad,det);
222 for(Int_t i=0;i<9;i++) rtd[i] = rt[i];
224 rm = new TRotMatrix(name,name,rtd);
226 cout << "EXCEPTION in new TRotMatrix" << endl;
230 //sprintf(name,"ND%1.1d2.2d2.2d",lay,lad,det);
232 nd = new TNode("SDD"," ",box,xg[0],xg[1],xg[2],rm);
234 cout << "EXCEPTION in new TNode" << endl;
237 nd->SetLineColor(kColorITS_SSD);
240 for(Int_t lay=5;lay<=5;lay++)
241 for(Int_t lad=1;lad<=gm->GetNladders(lay);lad++)
242 for(Int_t det=1;det<=gm->GetNdetectors(lay);det++){
244 box = new TBRIK ("ActiveSSD","Active volume of SSD","SSD SI DET",
247 cout << "EXCEPTION in box = new TBRIK" << endl;
250 gm->GetTrans(lay,lad,det,xg[0],xg[1],xg[2]);
251 gm->GetRotMatrix(lay,lad,det,rt);
252 //sprintf(name,"ROT%1.1d2.2d2.2d",lay,lad,det);
253 for(Int_t i=0;i<9;i++) rtd[i] = rt[i];
255 rm = new TRotMatrix(name,name,rtd);
257 cout << "EXCEPTION in new TRotMatrix" << endl;
261 //sprintf(name,"ND%1.1d2.2d2.2d",lay,lad,det);
263 nd = new TNode("SSD"," ",box,xg[0],xg[1],xg[2],rm);
265 cout << "EXCEPTION in new TNode" << endl;
268 nd->SetLineColor(kColorITS_SSD);
272 for(Int_t lay=6;lay<=6;lay++)
273 for(Int_t lad=1;lad<=gm->GetNladders(lay);lad++)
274 for(Int_t det=1;det<=gm->GetNdetectors(lay);det++){
276 box = new TBRIK ("ActiveSSD","Active volume of SSD","SSD SI DET",
279 cout << "EXCEPTION in box = new TBRIK" << endl;
283 gm->GetTrans(lay,lad,det,xg[0],xg[1],xg[2]);
284 gm->GetRotMatrix(lay,lad,det,rt);
285 //sprintf(name,"ROT%1.1d2.2d2.2d",lay,lad,det);
286 for(Int_t i=0;i<9;i++) rtd[i] = rt[i];
288 rm = new TRotMatrix(name,name,rtd);
290 cout << "EXCEPTION in new TRotMatrix" << endl;
294 //sprintf(name,"ND%1.1d2.2d2.2d",lay,lad,det);
296 nd = new TNode("SSD"," ",box,xg[0],xg[1],xg[2],rm);
298 cout << "EXCEPTION in new TNode" << endl;
301 nd->SetLineColor(kColorITS_SSD);
307 //_____________________________________________________________________________
308 void AliITSv5::CreateMaterials(){
309 ////////////////////////////////////////////////////////////////////////
310 // Read a file containing the materials for the ITS version 5.
311 ////////////////////////////////////////////////////////////////////////
314 filtmp = gSystem->ExpandPathName(fEuclidMaterial.Data());
316 FILE *file = fopen(filtmp,"r");
319 ReadEuclidMedia(filtmp);
321 Error("CreateMaterials"," THE MEDIA FILE %s DOES NOT EXIST !",filtmp);
325 //_____________________________________________________________________________
326 void AliITSv5::CreateGeometry(){
327 //////////////////////////////////////////////////////////////////////
328 // This is the geometry used for the ITS Pre-TDR and comes from an
329 // Euclid to Geant conversion. The only difference
330 // is in the details of the ITS supports. The detectors elements,
331 // detector numbering, and local and global reference frames are shown in
332 // the following figures.
335 <img src="picts/ITS/its1+2_convention_front_5.gif">
338 <font size=+2 color=red>
339 <p>This shows the front view of the SPDs.
342 <img src="picts/ITS/its1+2_convention_side_5.gif">
345 <font size=+2 color=red>
346 <p>This shows the perspective view of the SPDs.
348 <img src="picts/ITS/its1+2_tree.gif">
351 <font size=+2 color=red>
352 <p>This shows the geometry Tree for the SPDs.
357 <img src="picts/ITS/its3+4_convention_front_5.gif">
360 <font size=+2 color=red>
361 <p>This shows the front view of the SDDs.
364 <img src="picts/ITS/its3+4_convention_side_5.gif">
367 <font size=+2 color=red>
368 <p>This shows the perspective view of the SDDs.
370 <img src="picts/ITS/its3+4_tree.gif">
373 <font size=+2 color=red>
374 <p>This shows the geometry Tree for the SDDs.
379 <img src="picts/ITS/its5+6_convention_front_5.gif">
382 <font size=+2 color=red>
383 <p>This shows the front view of the SSDs.
386 <img src="picts/ITS/its5+6_convention_side_5.gif">
389 <font size=+2 color=red>
390 <p>This shows the perspective view of the SSDs.
393 <img src="picts/ITS/its5+6_tree.gif">
396 <font size=+2 color=red>
397 <p>This shows the geometry Tree for the SSDs.
402 <img src="picts/ITS/its_layer1-6_2.gif">
405 <font size=+2 color=red>
406 <p>This shows the front view of the whole ITS..
410 <img src="picts/ITS/its_layer1-6_1.gif">
413 <font size=+2 color=red>
414 <p>This shows the perspective view of the whole ITS..
418 <img src="picts/ITS/its1-6_tree.gif">
421 <font size=+2 color=red>
422 <p>This shows the geometry Tree for the whole ITS.
429 // Here are shown the details of the ITS support cones and services.
430 // First is a GEANT tree showing the organization of all of the volumes
431 // that make up the ITS supports and services.
434 <img src="picts/ITS/supports_tree.gif">
437 // What follows are a number of figures showing what these support
438 // structures look like.
442 <img src="picts/ITS/supports_3.gif">
445 <font size=+2 color=red>
446 <p>This shows the geometry of the supports for the Drift and Strip layers only.
450 <img src="picts/ITS/supports_2.gif">
453 <font size=+2 color=red>
454 <p>This shows the geometry of the supports for the Drift and Strip layers in front cut out.
458 <img src="picts/ITS/supports_1.gif">
461 <font size=+2 color=red>
462 <p>This shows the geometry of the supports for the Drift and Strip layers in a back cut out..
466 <img src="picts/ITS/suppssd.gif">
469 <font size=+2 color=red>
470 <p>This shows the geometry for the Strip layers supports.
474 <img src="picts/ITS/suppsdd.gif">
477 <font size=+2 color=red>
478 <p>This shows the geometry for the Drift layers supports.
484 // Read a file containing the geometry for the ITS version 5.
485 ////////////////////////////////////////////////////////////////////////
491 filtmp = gSystem->ExpandPathName(fEuclidGeometry.Data());
492 FILE *file = fopen(filtmp,"r");
496 printf("Ready to read Euclid geometry file\n");
497 ReadEuclid(fEuclidGeometry.Data(),topvol);
498 printf("Read in euclid geometries\n");
500 Error("CreateGeometry"," THE GEOM FILE %s DOES NOT EXIST !",
501 fEuclidGeometry.Data());
505 // Place the ITS ghost volume ITSV in its mother volume (ALIC) and make it
508 gMC->Gspos("ITSV",1,"ALIC",0,0,0,0,"ONLY");
510 // Outputs the geometry tree in the EUCLID/CAD format if requested to do so
513 gMC->WriteEuclid("ITSgeometry", "ITSV", 1, 5);
514 } // end if (fEuclidOut)
516 // read in the file containing the transformations for the active
517 // volumes for the ITS version 5. This is expected to be in a file
518 // ending in .det. This geometry is kept in the AliITSgeom class.
519 filtmp = gSystem->ExpandPathName(fEuclidGeometry.Data());
520 size = strlen(filtmp);
522 filtmp[size-3] = 'd'; // change from .euc to .det
523 filtmp[size-2] = 'e';
524 filtmp[size-1] = 't';
525 file = fopen(filtmp,"r");
526 if(file){ // if file exists use it to fill AliITSgeom structure.
528 printf("ready to read .det file %s\n",filtmp);
529 fITSgeom = new AliITSgeom(filtmp);
532 // fill AliITSgeom structure from geant structure just filled above
536 printf("finished with euclid geometrys\n");
538 //_____________________________________________________________________________
539 void AliITSv5::Init(){
540 ////////////////////////////////////////////////////////////////////////
541 // Initialise the ITS after it has been created.
542 ////////////////////////////////////////////////////////////////////////
546 fIdName = new char*[fId5N];
547 fIdSens = new Int_t[fId5N];
548 for(i=0;i<fId5N;i++) {
549 l = strlen(fId5Name[i]);
550 fIdName[i] = new char[l+1];
551 for(j=0;j<l;j++) fIdName[i][j] = fId5Name[i][j];
552 fIdName[i][l] = '\0'; // Null terminate this string.
559 //_____________________________________________________________________________
560 void AliITSv5::StepManager(){
561 ////////////////////////////////////////////////////////////////////////
562 // Called for every step in the ITS, then calles the AliITShit class
563 // creator with the information to be recoreded about that hit.
564 // The value of the macro ALIITSPRINTGEOM if set to 1 will allow the
565 // printing of information to a file which can be used to create a .det
566 // file read in by the routine CreateGeometry(). If set to 0 or any other
567 // value except 1, the default behavior, then no such file is created nor
568 // it the extra variables and the like used in the printing allocated.
569 ////////////////////////////////////////////////////////////////////////
574 TLorentzVector position, momentum;
575 TClonesArray &lhits = *fHits;
576 #if ALIITSPRINTGEOM==1
579 Float_t xl[3],xt[3],angl[6];
580 // Float_t par[20],att[20];
582 static Bool_t first=kTRUE,printit[6][50][50];
583 if(first){ for(copy1=0;copy1<6;copy1++)for(copy2=0;copy2<50;copy2++)
584 for(id=0;id<50;id++) printit[copy1][copy2][id] = kTRUE;
592 if(gMC->IsTrackInside()) vol[3] += 1;
593 if(gMC->IsTrackEntering()) vol[3] += 2;
594 if(gMC->IsTrackExiting()) vol[3] += 4;
595 if(gMC->IsTrackOut()) vol[3] += 8;
596 if(gMC->IsTrackDisappeared()) vol[3] += 16;
597 if(gMC->IsTrackStop()) vol[3] += 32;
598 if(gMC->IsTrackAlive()) vol[3] += 64;
600 // Fill hit structure.
601 if(!(gMC->TrackCharge())) return;
603 // Only entering charged tracks
604 if((id = gMC->CurrentVolID(copy)) == fIdSens[0]) {
606 id = gMC->CurrentVolOffID(0,copy);
607 //detector copy in the ladder = 1<->4 (ITS1)
609 gMC->CurrentVolOffID(1,copy1);
610 //ladder copy in the module = 1<->2 (I186)
611 gMC->CurrentVolOffID(2,copy2);
612 //module copy in the layer = 1<->10 (I132)
613 vol[2] = copy1+(copy2-1)*2;//# of ladders in one module = 2
614 } else if(id == fIdSens[1]){
616 id = gMC->CurrentVolOffID(0,copy);
617 //detector copy in the ladder = 1<->4 (ITS2)
619 gMC->CurrentVolOffID(1,copy1);
620 //ladder copy in the module = 1<->4 (I131)
621 gMC->CurrentVolOffID(2,copy2);
622 //module copy in the layer = 1<->10 (I132)
623 vol[2] = copy1+(copy2-1)*4;//# of ladders in one module = 4
624 } else if(id == fIdSens[2]){
626 id = gMC->CurrentVolOffID(1,copy);
627 //detector copy in the ladder = 1<->5 (ITS3 is inside I314)
629 id = gMC->CurrentVolOffID(2,copy);
630 //ladder copy in the layer = 1<->12 (I316)
632 } else if(id == fIdSens[3]){
634 id = gMC->CurrentVolOffID(1,copy);
635 //detector copy in the ladder = 1<->8 (ITS4 is inside I414)
637 id = gMC->CurrentVolOffID(2,copy);
638 //ladder copy in the layer = 1<->22 (I417)
640 }else if(id == fIdSens[4]){
642 id = gMC->CurrentVolOffID(1,copy);
643 //detector copy in the ladder = 1<->23 (ITS5 is inside I562)
645 id = gMC->CurrentVolOffID(2,copy);
646 //ladder copy in the layer = 1<->34 (I565)
648 }else if(id == fIdSens[5]){
650 id = gMC->CurrentVolOffID(1,copy);
651 //detector copy in the ladder = 1<->26 (ITS6 is inside I566)
653 id = gMC->CurrentVolOffID(2,copy);
654 //ladder copy in the layer = 1<->38 (I569)
657 return; // not an ITS volume?
658 } // end if/else if (gMC->CurentVolID(copy) == fIdSens[i])
660 gMC->TrackPosition(position);
661 gMC->TrackMomentum(momentum);
669 hits[7]=gMC->TrackTime();
670 // Fill hit structure with this new hit.
671 new(lhits[fNhits++]) AliITShit(fIshunt,gAlice->CurrentTrack(),vol,hits);
672 #if ALIITSPRINTGEOM==1
673 if(printit[vol[0]][vol[2]][vol[1]]){
674 printit[vol[0]][vol[2]][vol[1]] = kFALSE;
675 xl[0] = xl[1] = xl[2] = 0.0;
677 for(i=0;i<9;i++) mat[i] = 0.0;
678 mat[0] = mat[4] = mat[8] = 1.0; // default with identity matrix
681 gMC->Gdtom(xl,&(mat[0]),2);
684 gMC->Gdtom(xl,&(mat[3]),2);
687 gMC->Gdtom(xl,&(mat[6]),2);
689 angl[0] = TMath::ACos(mat[2]);
690 if(mat[2]==1.0) angl[0] = 0.0;
691 angl[1] = TMath::ATan2(mat[1],mat[0]);
692 if(angl[1]<0.0) angl[1] += 2.0*TMath::Pi();
694 angl[2] = TMath::ACos(mat[5]);
695 if(mat[5]==1.0) angl[2] = 0.0;
696 angl[3] = TMath::ATan2(mat[4],mat[3]);
697 if(angl[3]<0.0) angl[3] += 2.0*TMath::Pi();
699 angl[4] = TMath::ACos(mat[8]);
700 if(mat[8]==1.0) angl[4] = 0.0;
701 angl[5] = TMath::ATan2(mat[7],mat[6]);
702 if(angl[5]<0.0) angl[5] += 2.0*TMath::Pi();
704 for(i=0;i<6;i++) angl[i] *= 180.0/TMath::Pi(); // degrees
705 fp = fopen("ITSgeometry_v5.det","a");
706 fprintf(fp,"%2d %2d %2d %9e %9e %9e %9e %9e %9e %9e %9e %9e ",
707 vol[0],vol[2],vol[1], // layer ladder detector
708 xt[0],xt[1],xt[2], // Translation vector
709 angl[0],angl[1],angl[2],angl[3],angl[4],angl[5] // Geant rotaion
712 fprintf(fp,"%9e %9e %9e %9e %9e %9e %9e %9e %9e",
713 mat[0],mat[1],mat[2],mat[3],mat[4],mat[5],mat[6],mat[7],mat[8]
714 ); // Adding the rotation matrix.
717 } // end if printit[layer][ladder][detector]
721 //____________________________________________________________________________
722 void AliITSv5::Streamer(TBuffer &R__b){
723 ////////////////////////////////////////////////////////////////////////
724 // A dummy Streamer function for this class AliITSv5. By default it
725 // only streams the AliITS class as it is required. Since this class
726 // dosen't contain any "real" data to be saved, it doesn't.
727 ////////////////////////////////////////////////////////////////////////
728 if (R__b.IsReading()) {
729 Version_t R__v = R__b.ReadVersion();
731 AliITS::Streamer(R__b);
732 // This information does not need to be read. It is "hard wired"
733 // into this class via its creators.
735 //R__b.ReadArray(fId5Name);
739 R__b.WriteVersion(AliITSv5::IsA());
740 AliITS::Streamer(R__b);
741 // This information does not need to be saved. It is "hard wired"
742 // into this class via its creators.
744 //R__b.WriteArray(fId5Name, __COUNTER__);
745 } // end if R__b.IsReading()