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.17 2000/05/10 19:57:26 nilsen
19 Fixed problem with display.C and ITS versions v1 and v3.
21 Revision 1.16 2000/03/05 00:11:03 nilsen
24 Revision 1.15 2000/02/23 16:25:21 fca
25 AliVMC and AliGeant3 classes introduced
26 ReadEuclid moved from AliRun to AliModule
28 Revision 1.14.4.3 2000/03/04 23:46:38 nilsen
29 Fixed up the comments/documentation.
31 Revision 1.14.4.2 2000/03/02 21:53:02 nilsen
32 To make it compatable with the changes in AliRun/AliModule.
34 Revision 1.14.4.1 2000/01/12 19:03:33 nilsen
35 This is the version of the files after the merging done in December 1999.
36 See the ReadMe110100.txt file for details
38 Revision 1.14 1999/10/22 08:16:49 fca
39 Correct destructors, thanks to I.Hrivnacova
41 Revision 1.13 1999/10/06 10:15:19 fca
42 Correct bug in allocation of layer name and add destructor
44 Revision 1.12 1999/10/05 08:05:09 fca
45 Minor corrections for uninitialised variables.
47 Revision 1.11 1999/09/29 09:24:20 fca
48 Introduction of the Copyright and cvs Log
52 ///////////////////////////////////////////////////////////////////////////////
54 // Inner Traking System version 5
55 // This class contains the base procedures for the Inner Tracking System
57 // Authors: R. Barbera, B. S. Nilsen.
59 // Created September 17 1999.
61 ///////////////////////////////////////////////////////////////////////////////
63 // See AliITSv5::StepManager().
64 #define ALIITSPRINTGEOM 0 // default. don't print out gemetry information
65 //#define ALIITSPRINTGEOM 1 // print out geometry information
70 #include <TGeometry.h>
76 #if ALIITSPRINTGEOM==1
77 #include "../TGeant3/TGeant3.h"
79 #include "AliITShit.h"
82 #include "AliITSgeom.h"
86 //_____________________________________________________________________________
87 AliITSv5::AliITSv5() {
88 ////////////////////////////////////////////////////////////////////////
89 // Standard default constructor for the ITS version 5.
90 ////////////////////////////////////////////////////////////////////////
92 fId5Name = new char*[fId5N];
100 //_____________________________________________________________________________
101 AliITSv5::~AliITSv5() {
102 ////////////////////////////////////////////////////////////////////////
103 // Standard destructor for the ITS version 5.
104 ////////////////////////////////////////////////////////////////////////
107 //_____________________________________________________________________________
108 AliITSv5::AliITSv5(const char *name, const char *title) : AliITS(name, title){
109 ////////////////////////////////////////////////////////////////////////
110 // Standard constructor for the ITS version 5.
111 ////////////////////////////////////////////////////////////////////////
113 fId5Name = new char*[fId5N];
114 fId5Name[0] = "ITS1";
115 fId5Name[1] = "ITS2";
116 fId5Name[2] = "ITS3";
117 fId5Name[3] = "ITS4";
118 fId5Name[4] = "ITS5";
119 fId5Name[5] = "ITS6";
121 fEuclidMaterial = "$ALICE_ROOT/Euclid/ITSgeometry_5.tme";
122 fEuclidGeometry = "$ALICE_ROOT/Euclid/ITSgeometry_5.euc";
124 //_____________________________________________________________________________
125 void AliITSv5::BuildGeometry(){
127 // Build ITS TNODE geometry for event display using detailed geometry.
128 // This function builds a simple ITS geometry used by the ROOT macro
133 //const int kColorITS_SPD=kRed;
134 //const int kColorITS_SDD=kGreen;
135 const int kColorITS_SSD=kBlue;
137 Top=gAlice->GetGeometry()->GetNode("alice");
138 AliITSgeom *gm = this->GetITSgeom();
147 //TCanvas *c1 = new TCanvas("c1","ITS");
149 for(lay=1;lay<=2;lay++)
150 for(lad=1;lad<=gm->GetNladders(lay);lad++)
151 for(det=1;det<=gm->GetNdetectors(lay);det++){
153 box = new TBRIK ("ActiveSPD","Active volume of SPD","SPD SI DET",
156 cout << "EXCEPTION in box = new TBRIK" << endl;
159 gm->GetTrans(lay,lad,det,xg[0],xg[1],xg[2]);
160 gm->GetRotMatrix(lay,lad,det,rt);
161 //sprintf(name,"ROT%1.1d2.2d2.2d",lay,lad,det);
162 for(Int_t i=0;i<9;i++) rtd[i] = rt[i];
164 rm = new TRotMatrix(name,name,rtd);
166 cout << "EXCEPTION in new TRotMatrix" << endl;
170 //sprintf(name,"ND%1.1d2.2d2.2d",lay,lad,det);
172 nd = new TNode("SPD"," ",box,xg[0],xg[1],xg[2],rm);
174 cout << "EXCEPTION in new TNode" << endl;
177 nd->SetLineColor(kColorITS_SSD);
181 for(lay=3;lay<=3;lay++)
182 for(lad=1;lad<=gm->GetNladders(lay);lad++)
183 for(det=1;det<=gm->GetNdetectors(lay);det++){
185 box = new TBRIK ("ActiveSDD","Active volume of SDD","SDD SI DET",
188 cout << "EXCEPTION in box = new TBRIK" << endl;
191 gm->GetTrans(lay,lad,det,xg[0],xg[1],xg[2]);
192 gm->GetRotMatrix(lay,lad,det,rt);
193 //sprintf(name,"ROT%1.1d2.2d2.2d",lay,lad,det);
194 for(Int_t i=0;i<9;i++) rtd[i] = rt[i];
196 rm = new TRotMatrix(name,name,rtd);
198 cout << "EXCEPTION in new TRotMatrix" << endl;
202 //sprintf(name,"ND%1.1d2.2d2.2d",lay,lad,det);
204 nd = new TNode("SDD"," ",box,xg[0],xg[1],xg[2],rm);
206 cout << "EXCEPTION in new TNode" << endl;
209 nd->SetLineColor(kColorITS_SSD);
213 for(lay=4;lay<=4;lay++)
214 for(lad=1;lad<=gm->GetNladders(lay);lad++)
215 for(det=1;det<=gm->GetNdetectors(lay);det++){
217 box = new TBRIK ("ActiveSDD","Active volume of SDD","SDD SI DET",
220 cout << "EXCEPTION in box = new TBRIK" << endl;
223 gm->GetTrans(lay,lad,det,xg[0],xg[1],xg[2]);
224 gm->GetRotMatrix(lay,lad,det,rt);
225 //sprintf(name,"ROT%1.1d2.2d2.2d",lay,lad,det);
226 for(Int_t i=0;i<9;i++) rtd[i] = rt[i];
228 rm = new TRotMatrix(name,name,rtd);
230 cout << "EXCEPTION in new TRotMatrix" << endl;
234 //sprintf(name,"ND%1.1d2.2d2.2d",lay,lad,det);
236 nd = new TNode("SDD"," ",box,xg[0],xg[1],xg[2],rm);
238 cout << "EXCEPTION in new TNode" << endl;
241 nd->SetLineColor(kColorITS_SSD);
244 for(lay=5;lay<=5;lay++)
245 for(lad=1;lad<=gm->GetNladders(lay);lad++)
246 for(det=1;det<=gm->GetNdetectors(lay);det++){
248 box = new TBRIK ("ActiveSSD","Active volume of SSD","SSD SI DET",
251 cout << "EXCEPTION in box = new TBRIK" << endl;
254 gm->GetTrans(lay,lad,det,xg[0],xg[1],xg[2]);
255 gm->GetRotMatrix(lay,lad,det,rt);
256 //sprintf(name,"ROT%1.1d2.2d2.2d",lay,lad,det);
257 for(Int_t i=0;i<9;i++) rtd[i] = rt[i];
259 rm = new TRotMatrix(name,name,rtd);
261 cout << "EXCEPTION in new TRotMatrix" << endl;
265 //sprintf(name,"ND%1.1d2.2d2.2d",lay,lad,det);
267 nd = new TNode("SSD"," ",box,xg[0],xg[1],xg[2],rm);
269 cout << "EXCEPTION in new TNode" << endl;
272 nd->SetLineColor(kColorITS_SSD);
276 for(lay=6;lay<=6;lay++)
277 for(lad=1;lad<=gm->GetNladders(lay);lad++)
278 for(det=1;det<=gm->GetNdetectors(lay);det++){
280 box = new TBRIK ("ActiveSSD","Active volume of SSD","SSD SI DET",
283 cout << "EXCEPTION in box = new TBRIK" << endl;
287 gm->GetTrans(lay,lad,det,xg[0],xg[1],xg[2]);
288 gm->GetRotMatrix(lay,lad,det,rt);
289 //sprintf(name,"ROT%1.1d2.2d2.2d",lay,lad,det);
290 for(i=0;i<9;i++) rtd[i] = rt[i];
292 rm = new TRotMatrix(name,name,rtd);
294 cout << "EXCEPTION in new TRotMatrix" << endl;
298 //sprintf(name,"ND%1.1d2.2d2.2d",lay,lad,det);
300 nd = new TNode("SSD"," ",box,xg[0],xg[1],xg[2],rm);
302 cout << "EXCEPTION in new TNode" << endl;
305 nd->SetLineColor(kColorITS_SSD);
311 //_____________________________________________________________________________
312 void AliITSv5::CreateMaterials(){
313 ////////////////////////////////////////////////////////////////////////
314 // Read a file containing the materials for the ITS version 5.
315 ////////////////////////////////////////////////////////////////////////
318 filtmp = gSystem->ExpandPathName(fEuclidMaterial.Data());
320 FILE *file = fopen(filtmp,"r");
323 ReadEuclidMedia(filtmp);
325 Error("CreateMaterials"," THE MEDIA FILE %s DOES NOT EXIST !",filtmp);
329 //_____________________________________________________________________________
330 void AliITSv5::CreateGeometry(){
331 //////////////////////////////////////////////////////////////////////
332 // This is the geometry used for the ITS Pre-TDR and comes from an
333 // Euclid to Geant conversion. The only difference
334 // is in the details of the ITS supports. The detectors elements,
335 // detector numbering, and local and global reference frames are shown in
336 // the following figures.
339 <img src="picts/ITS/its1+2_convention_front_5.gif">
342 <font size=+2 color=red>
343 <p>This shows the front view of the SPDs.
346 <img src="picts/ITS/its1+2_convention_side_5.gif">
349 <font size=+2 color=red>
350 <p>This shows the perspective view of the SPDs.
352 <img src="picts/ITS/its1+2_tree.gif">
355 <font size=+2 color=red>
356 <p>This shows the geometry Tree for the SPDs.
361 <img src="picts/ITS/its3+4_convention_front_5.gif">
364 <font size=+2 color=red>
365 <p>This shows the front view of the SDDs.
368 <img src="picts/ITS/its3+4_convention_side_5.gif">
371 <font size=+2 color=red>
372 <p>This shows the perspective view of the SDDs.
374 <img src="picts/ITS/its3+4_tree.gif">
377 <font size=+2 color=red>
378 <p>This shows the geometry Tree for the SDDs.
383 <img src="picts/ITS/its5+6_convention_front_5.gif">
386 <font size=+2 color=red>
387 <p>This shows the front view of the SSDs.
390 <img src="picts/ITS/its5+6_convention_side_5.gif">
393 <font size=+2 color=red>
394 <p>This shows the perspective view of the SSDs.
397 <img src="picts/ITS/its5+6_tree.gif">
400 <font size=+2 color=red>
401 <p>This shows the geometry Tree for the SSDs.
406 <img src="picts/ITS/its_layer1-6_2.gif">
409 <font size=+2 color=red>
410 <p>This shows the front view of the whole ITS..
414 <img src="picts/ITS/its_layer1-6_1.gif">
417 <font size=+2 color=red>
418 <p>This shows the perspective view of the whole ITS..
422 <img src="picts/ITS/its1-6_tree.gif">
425 <font size=+2 color=red>
426 <p>This shows the geometry Tree for the whole ITS.
433 // Here are shown the details of the ITS support cones and services.
434 // First is a GEANT tree showing the organization of all of the volumes
435 // that make up the ITS supports and services.
438 <img src="picts/ITS/supports_tree.gif">
441 // What follows are a number of figures showing what these support
442 // structures look like.
446 <img src="picts/ITS/supports_3.gif">
449 <font size=+2 color=red>
450 <p>This shows the geometry of the supports for the Drift and Strip layers only.
454 <img src="picts/ITS/supports_2.gif">
457 <font size=+2 color=red>
458 <p>This shows the geometry of the supports for the Drift and Strip layers in front cut out.
462 <img src="picts/ITS/supports_1.gif">
465 <font size=+2 color=red>
466 <p>This shows the geometry of the supports for the Drift and Strip layers in a back cut out..
470 <img src="picts/ITS/suppssd.gif">
473 <font size=+2 color=red>
474 <p>This shows the geometry for the Strip layers supports.
478 <img src="picts/ITS/suppsdd.gif">
481 <font size=+2 color=red>
482 <p>This shows the geometry for the Drift layers supports.
488 // Read a file containing the geometry for the ITS version 5.
489 ////////////////////////////////////////////////////////////////////////
495 filtmp = gSystem->ExpandPathName(fEuclidGeometry.Data());
496 FILE *file = fopen(filtmp,"r");
500 printf("Ready to read Euclid geometry file\n");
501 ReadEuclid(fEuclidGeometry.Data(),topvol);
502 printf("Read in euclid geometries\n");
504 Error("CreateGeometry"," THE GEOM FILE %s DOES NOT EXIST !",
505 fEuclidGeometry.Data());
509 // Place the ITS ghost volume ITSV in its mother volume (ALIC) and make it
512 gMC->Gspos("ITSV",1,"ALIC",0,0,0,0,"ONLY");
514 // Outputs the geometry tree in the EUCLID/CAD format if requested to do so
517 gMC->WriteEuclid("ITSgeometry", "ITSV", 1, 5);
518 } // end if (fEuclidOut)
520 // read in the file containing the transformations for the active
521 // volumes for the ITS version 5. This is expected to be in a file
522 // ending in .det. This geometry is kept in the AliITSgeom class.
523 filtmp = gSystem->ExpandPathName(fEuclidGeometry.Data());
524 size = strlen(filtmp);
526 filtmp[size-3] = 'd'; // change from .euc to .det
527 filtmp[size-2] = 'e';
528 filtmp[size-1] = 't';
529 file = fopen(filtmp,"r");
530 if(file){ // if file exists use it to fill AliITSgeom structure.
532 printf("ready to read .det file %s\n",filtmp);
533 fITSgeom = new AliITSgeom(filtmp);
536 // fill AliITSgeom structure from geant structure just filled above
540 printf("finished with euclid geometrys\n");
542 //_____________________________________________________________________________
543 void AliITSv5::Init(){
544 ////////////////////////////////////////////////////////////////////////
545 // Initialise the ITS after it has been created.
546 ////////////////////////////////////////////////////////////////////////
550 fIdName = new char*[fId5N];
551 fIdSens = new Int_t[fId5N];
552 for(i=0;i<fId5N;i++) {
553 l = strlen(fId5Name[i]);
554 fIdName[i] = new char[l+1];
555 for(j=0;j<l;j++) fIdName[i][j] = fId5Name[i][j];
556 fIdName[i][l] = '\0'; // Null terminate this string.
563 //_____________________________________________________________________________
564 void AliITSv5::StepManager(){
565 ////////////////////////////////////////////////////////////////////////
566 // Called for every step in the ITS, then calles the AliITShit class
567 // creator with the information to be recoreded about that hit.
568 // The value of the macro ALIITSPRINTGEOM if set to 1 will allow the
569 // printing of information to a file which can be used to create a .det
570 // file read in by the routine CreateGeometry(). If set to 0 or any other
571 // value except 1, the default behavior, then no such file is created nor
572 // it the extra variables and the like used in the printing allocated.
573 ////////////////////////////////////////////////////////////////////////
578 TLorentzVector position, momentum;
579 TClonesArray &lhits = *fHits;
580 #if ALIITSPRINTGEOM==1
583 Float_t xl[3],xt[3],angl[6];
584 // Float_t par[20],att[20];
586 static Bool_t first=kTRUE,printit[6][50][50];
587 if(first){ for(copy1=0;copy1<6;copy1++)for(copy2=0;copy2<50;copy2++)
588 for(id=0;id<50;id++) printit[copy1][copy2][id] = kTRUE;
596 if(gMC->IsTrackInside()) vol[3] += 1;
597 if(gMC->IsTrackEntering()) vol[3] += 2;
598 if(gMC->IsTrackExiting()) vol[3] += 4;
599 if(gMC->IsTrackOut()) vol[3] += 8;
600 if(gMC->IsTrackDisappeared()) vol[3] += 16;
601 if(gMC->IsTrackStop()) vol[3] += 32;
602 if(gMC->IsTrackAlive()) vol[3] += 64;
604 // Fill hit structure.
605 if(!(gMC->TrackCharge())) return;
607 // Only entering charged tracks
608 if((id = gMC->CurrentVolID(copy)) == fIdSens[0]) {
610 id = gMC->CurrentVolOffID(0,copy);
611 //detector copy in the ladder = 1<->4 (ITS1)
613 gMC->CurrentVolOffID(1,copy1);
614 //ladder copy in the module = 1<->2 (I186)
615 gMC->CurrentVolOffID(2,copy2);
616 //module copy in the layer = 1<->10 (I132)
617 vol[2] = copy1+(copy2-1)*2;//# of ladders in one module = 2
618 } else if(id == fIdSens[1]){
620 id = gMC->CurrentVolOffID(0,copy);
621 //detector copy in the ladder = 1<->4 (ITS2)
623 gMC->CurrentVolOffID(1,copy1);
624 //ladder copy in the module = 1<->4 (I131)
625 gMC->CurrentVolOffID(2,copy2);
626 //module copy in the layer = 1<->10 (I132)
627 vol[2] = copy1+(copy2-1)*4;//# of ladders in one module = 4
628 } else if(id == fIdSens[2]){
630 id = gMC->CurrentVolOffID(1,copy);
631 //detector copy in the ladder = 1<->5 (ITS3 is inside I314)
633 id = gMC->CurrentVolOffID(2,copy);
634 //ladder copy in the layer = 1<->12 (I316)
636 } else if(id == fIdSens[3]){
638 id = gMC->CurrentVolOffID(1,copy);
639 //detector copy in the ladder = 1<->8 (ITS4 is inside I414)
641 id = gMC->CurrentVolOffID(2,copy);
642 //ladder copy in the layer = 1<->22 (I417)
644 }else if(id == fIdSens[4]){
646 id = gMC->CurrentVolOffID(1,copy);
647 //detector copy in the ladder = 1<->23 (ITS5 is inside I562)
649 id = gMC->CurrentVolOffID(2,copy);
650 //ladder copy in the layer = 1<->34 (I565)
652 }else if(id == fIdSens[5]){
654 id = gMC->CurrentVolOffID(1,copy);
655 //detector copy in the ladder = 1<->26 (ITS6 is inside I566)
657 id = gMC->CurrentVolOffID(2,copy);
658 //ladder copy in the layer = 1<->38 (I569)
661 return; // not an ITS volume?
662 } // end if/else if (gMC->CurentVolID(copy) == fIdSens[i])
664 gMC->TrackPosition(position);
665 gMC->TrackMomentum(momentum);
673 hits[7]=gMC->TrackTime();
674 // Fill hit structure with this new hit.
675 new(lhits[fNhits++]) AliITShit(fIshunt,gAlice->CurrentTrack(),vol,hits);
676 #if ALIITSPRINTGEOM==1
677 if(printit[vol[0]][vol[2]][vol[1]]){
678 printit[vol[0]][vol[2]][vol[1]] = kFALSE;
679 xl[0] = xl[1] = xl[2] = 0.0;
681 for(i=0;i<9;i++) mat[i] = 0.0;
682 mat[0] = mat[4] = mat[8] = 1.0; // default with identity matrix
685 gMC->Gdtom(xl,&(mat[0]),2);
688 gMC->Gdtom(xl,&(mat[3]),2);
691 gMC->Gdtom(xl,&(mat[6]),2);
693 angl[0] = TMath::ACos(mat[2]);
694 if(mat[2]==1.0) angl[0] = 0.0;
695 angl[1] = TMath::ATan2(mat[1],mat[0]);
696 if(angl[1]<0.0) angl[1] += 2.0*TMath::Pi();
698 angl[2] = TMath::ACos(mat[5]);
699 if(mat[5]==1.0) angl[2] = 0.0;
700 angl[3] = TMath::ATan2(mat[4],mat[3]);
701 if(angl[3]<0.0) angl[3] += 2.0*TMath::Pi();
703 angl[4] = TMath::ACos(mat[8]);
704 if(mat[8]==1.0) angl[4] = 0.0;
705 angl[5] = TMath::ATan2(mat[7],mat[6]);
706 if(angl[5]<0.0) angl[5] += 2.0*TMath::Pi();
708 for(i=0;i<6;i++) angl[i] *= 180.0/TMath::Pi(); // degrees
709 fp = fopen("ITSgeometry_v5.det","a");
710 fprintf(fp,"%2d %2d %2d %9e %9e %9e %9e %9e %9e %9e %9e %9e ",
711 vol[0],vol[2],vol[1], // layer ladder detector
712 xt[0],xt[1],xt[2], // Translation vector
713 angl[0],angl[1],angl[2],angl[3],angl[4],angl[5] // Geant rotaion
716 fprintf(fp,"%9e %9e %9e %9e %9e %9e %9e %9e %9e",
717 mat[0],mat[1],mat[2],mat[3],mat[4],mat[5],mat[6],mat[7],mat[8]
718 ); // Adding the rotation matrix.
721 } // end if printit[layer][ladder][detector]
725 //____________________________________________________________________________
726 void AliITSv5::Streamer(TBuffer &R__b){
727 ////////////////////////////////////////////////////////////////////////
728 // A dummy Streamer function for this class AliITSv5. By default it
729 // only streams the AliITS class as it is required. Since this class
730 // dosen't contain any "real" data to be saved, it doesn't.
731 ////////////////////////////////////////////////////////////////////////
732 if (R__b.IsReading()) {
733 Version_t R__v = R__b.ReadVersion();
735 AliITS::Streamer(R__b);
736 // This information does not need to be read. It is "hard wired"
737 // into this class via its creators.
739 //R__b.ReadArray(fId5Name);
743 R__b.WriteVersion(AliITSv5::IsA());
744 AliITS::Streamer(R__b);
745 // This information does not need to be saved. It is "hard wired"
746 // into this class via its creators.
748 //R__b.WriteArray(fId5Name, __COUNTER__);
749 } // end if R__b.IsReading()