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.25 2000/10/05 20:50:00 nilsen
19 Now using root generated streamers.
21 Revision 1.14.4.12 2000/10/02 16:04:03 barbera
22 Forward declarations added
24 Revision 1.22 2000/07/10 16:07:19 fca
25 Release version of ITS code
27 Revision 1.14.4.4 2000/05/19 10:10:21 nilsen
28 fix for bug with HP and Sun unix + fix for event display in ITS-working branch
30 Revision 1.14.4.3 2000/03/04 23:46:38 nilsen
31 Fixed up the comments/documentation.
33 Revision 1.14.4.2 2000/03/02 21:53:02 nilsen
34 To make it compatable with the changes in AliRun/AliModule.
36 Revision 1.14.4.1 2000/01/12 19:03:33 nilsen
37 This is the version of the files after the merging done in December 1999.
38 See the ReadMe110100.txt file for details
40 Revision 1.14 1999/10/22 08:16:49 fca
41 Correct destructors, thanks to I.Hrivnacova
43 Revision 1.13 1999/10/06 10:15:19 fca
44 Correct bug in allocation of layer name and add destructor
46 Revision 1.12 1999/10/05 08:05:09 fca
47 Minor corrections for uninitialised variables.
49 Revision 1.11 1999/09/29 09:24:20 fca
50 Introduction of the Copyright and cvs Log
54 ///////////////////////////////////////////////////////////////////////////////
56 // Inner Traking System version 5 with symmetric services
57 // This class contains the base procedures for the Inner Tracking System
59 // Authors: R. Barbera, B. S. Nilsen.
61 // Created October 7 2000.
63 ///////////////////////////////////////////////////////////////////////////////
65 // See AliITSv5symm::StepManager().
66 #define ALIITSPRINTGEOM 0 // default. don't print out gemetry information
67 //#define ALIITSPRINTGEOM 1 // print out geometry information
72 #include <TGeometry.h>
75 #include <TFile.h> // only required for Tracking function?
77 #include <TObjArray.h>
78 #include <TObjString.h>
79 #include <TClonesArray.h>
85 #if ALIITSPRINTGEOM==1
86 #include "../TGeant3/TGeant3.h"
88 #include "AliITShit.h"
90 #include "AliITSv5symm.h"
91 #include "AliITSgeom.h"
93 ClassImp(AliITSv5symm)
95 //_____________________________________________________________________________
96 AliITSv5symm::AliITSv5symm() {
97 ////////////////////////////////////////////////////////////////////////
98 // Standard default constructor for the ITS version 5.
99 ////////////////////////////////////////////////////////////////////////
107 //____________________________________________________________________________
108 AliITSv5symm::AliITSv5symm(const AliITSv5symm &source){
109 ////////////////////////////////////////////////////////////////////////
110 // Copy Constructor for ITS version 5.
111 ////////////////////////////////////////////////////////////////////////
112 if(&source == this) return;
113 printf("Not allowed to copy AliITSv5symm\n");
116 //_____________________________________________________________________________
117 AliITSv5symm& AliITSv5symm::operator=(const AliITSv5symm &source){
118 ////////////////////////////////////////////////////////////////////////
119 // Assignment operator for the ITS version 1.
120 ////////////////////////////////////////////////////////////////////////
122 if(&source == this) return *this;
123 printf("Not allowed to copy AliITSv3\n");
127 //_____________________________________________________________________________
128 AliITSv5symm::~AliITSv5symm() {
129 ////////////////////////////////////////////////////////////////////////
130 // Standard destructor for the ITS version 5.
131 ////////////////////////////////////////////////////////////////////////
133 //_____________________________________________________________________________
134 AliITSv5symm::AliITSv5symm(const char *name, const char *title) : AliITS(name, title){
135 /////////////////////////////////////////////////////////////////////////////
136 // Standard constructor for the ITS version 5 with symmetrical services.
137 /////////////////////////////////////////////////////////////////////////////
140 // TObjArray of TObjStrings
141 fIdName = new TObjArray(fIdN);
142 fIdName->AddAt(new TObjString("ITS1"),0);
143 fIdName->AddAt(new TObjString("ITS2"),1);
144 fIdName->AddAt(new TObjString("ITS3"),2);
145 fIdName->AddAt(new TObjString("ITS4"),3);
146 fIdName->AddAt(new TObjString("ITS5"),4);
147 fIdName->AddAt(new TObjString("ITS6"),5);
149 // Array of TStrings.
150 fIdName = new TString[fIdN];
157 fIdSens = new Int_t[fIdN];
158 for (Int_t i=0;i<fIdN;i++) fIdSens[i] = 0;
162 fEuclidMaterial = "$ALICE_ROOT/Euclid/ITSgeometry_5symm.tme";
163 fEuclidGeometry = "$ALICE_ROOT/Euclid/ITSgeometry_5symm.euc";
165 //_____________________________________________________________________________
166 void AliITSv5symm::BuildGeometry(){
167 ////////////////////////////////////////////////////////////////////////
168 // Geometry builder for the ITS version 5.
169 ////////////////////////////////////////////////////////////////////////
171 // Build ITS TNODE geometry for event display using detailed geometry.
172 // This function builds a simple ITS geometry used by the ROOT macro
177 //const int kColorITSSPD=kRed;
178 //const int kColorITSSDD=kGreen;
179 const int kColorITSSSD=kBlue;
181 top=gAlice->GetGeometry()->GetNode("alice");
182 AliITSgeom *gm = this->GetITSgeom();
191 //TCanvas *c1 = new TCanvas("c1","ITS");
193 for(lay=1;lay<=2;lay++)
194 for(lad=1;lad<=gm->GetNladders(lay);lad++)
195 for(det=1;det<=gm->GetNdetectors(lay);det++){
197 box = new TBRIK ("ActiveSPD","Active volume of SPD","SPD SI DET",
200 cout << "EXCEPTION in box = new TBRIK" << endl;
203 gm->GetTrans(lay,lad,det,xg[0],xg[1],xg[2]);
204 gm->GetRotMatrix(lay,lad,det,rt);
205 //sprintf(name,"ROT%1.1d2.2d2.2d",lay,lad,det);
206 for(i=0;i<9;i++) rtd[i] = rt[i];
208 rm = new TRotMatrix(name,name,rtd);
210 cout << "EXCEPTION in new TRotMatrix" << endl;
214 //sprintf(name,"ND%1.1d2.2d2.2d",lay,lad,det);
216 nd = new TNode("SPD"," ",box,xg[0],xg[1],xg[2],rm);
218 cout << "EXCEPTION in new TNode" << endl;
221 nd->SetLineColor(kColorITSSSD);
225 for(lay=3;lay<=3;lay++)
226 for(lad=1;lad<=gm->GetNladders(lay);lad++)
227 for(det=1;det<=gm->GetNdetectors(lay);det++){
229 box = new TBRIK ("ActiveSDD","Active volume of SDD","SDD SI DET",
232 cout << "EXCEPTION in box = new TBRIK" << endl;
235 gm->GetTrans(lay,lad,det,xg[0],xg[1],xg[2]);
236 gm->GetRotMatrix(lay,lad,det,rt);
237 //sprintf(name,"ROT%1.1d2.2d2.2d",lay,lad,det);
238 for(i=0;i<9;i++) rtd[i] = rt[i];
240 rm = new TRotMatrix(name,name,rtd);
242 cout << "EXCEPTION in new TRotMatrix" << endl;
246 //sprintf(name,"ND%1.1d2.2d2.2d",lay,lad,det);
248 nd = new TNode("SDD"," ",box,xg[0],xg[1],xg[2],rm);
250 cout << "EXCEPTION in new TNode" << endl;
253 nd->SetLineColor(kColorITSSSD);
257 for(lay=4;lay<=4;lay++)
258 for(lad=1;lad<=gm->GetNladders(lay);lad++)
259 for(det=1;det<=gm->GetNdetectors(lay);det++){
261 box = new TBRIK ("ActiveSDD","Active volume of SDD","SDD SI DET",
264 cout << "EXCEPTION in box = new TBRIK" << endl;
267 gm->GetTrans(lay,lad,det,xg[0],xg[1],xg[2]);
268 gm->GetRotMatrix(lay,lad,det,rt);
269 //sprintf(name,"ROT%1.1d2.2d2.2d",lay,lad,det);
270 for(i=0;i<9;i++) rtd[i] = rt[i];
272 rm = new TRotMatrix(name,name,rtd);
274 cout << "EXCEPTION in new TRotMatrix" << endl;
278 //sprintf(name,"ND%1.1d2.2d2.2d",lay,lad,det);
280 nd = new TNode("SDD"," ",box,xg[0],xg[1],xg[2],rm);
282 cout << "EXCEPTION in new TNode" << endl;
285 nd->SetLineColor(kColorITSSSD);
288 for(lay=5;lay<=5;lay++)
289 for(lad=1;lad<=gm->GetNladders(lay);lad++)
290 for(det=1;det<=gm->GetNdetectors(lay);det++){
292 box = new TBRIK ("ActiveSSD","Active volume of SSD","SSD SI DET",
295 cout << "EXCEPTION in box = new TBRIK" << endl;
298 gm->GetTrans(lay,lad,det,xg[0],xg[1],xg[2]);
299 gm->GetRotMatrix(lay,lad,det,rt);
300 //sprintf(name,"ROT%1.1d2.2d2.2d",lay,lad,det);
301 for(i=0;i<9;i++) rtd[i] = rt[i];
303 rm = new TRotMatrix(name,name,rtd);
305 cout << "EXCEPTION in new TRotMatrix" << endl;
309 //sprintf(name,"ND%1.1d2.2d2.2d",lay,lad,det);
311 nd = new TNode("SSD"," ",box,xg[0],xg[1],xg[2],rm);
313 cout << "EXCEPTION in new TNode" << endl;
316 nd->SetLineColor(kColorITSSSD);
320 for(lay=6;lay<=6;lay++)
321 for(lad=1;lad<=gm->GetNladders(lay);lad++)
322 for(det=1;det<=gm->GetNdetectors(lay);det++){
324 box = new TBRIK ("ActiveSSD","Active volume of SSD","SSD SI DET",
327 cout << "EXCEPTION in box = new TBRIK" << endl;
331 gm->GetTrans(lay,lad,det,xg[0],xg[1],xg[2]);
332 gm->GetRotMatrix(lay,lad,det,rt);
333 //sprintf(name,"ROT%1.1d2.2d2.2d",lay,lad,det);
334 for(i=0;i<9;i++) rtd[i] = rt[i];
336 rm = new TRotMatrix(name,name,rtd);
338 cout << "EXCEPTION in new TRotMatrix" << endl;
342 //sprintf(name,"ND%1.1d2.2d2.2d",lay,lad,det);
344 nd = new TNode("SSD"," ",box,xg[0],xg[1],xg[2],rm);
346 cout << "EXCEPTION in new TNode" << endl;
349 nd->SetLineColor(kColorITSSSD);
355 //_____________________________________________________________________________
356 void AliITSv5symm::CreateMaterials(){
357 ////////////////////////////////////////////////////////////////////////
358 // Read a file containing the materials for the ITS version 5.
359 ////////////////////////////////////////////////////////////////////////
362 filtmp = gSystem->ExpandPathName(fEuclidMaterial.Data());
364 FILE *file = fopen(filtmp,"r");
367 ReadEuclidMedia(filtmp);
369 Error("CreateMaterials"," THE MEDIA FILE %s DOES NOT EXIST !",filtmp);
373 //_____________________________________________________________________________
374 void AliITSv5symm::CreateGeometry(){
375 //////////////////////////////////////////////////////////////////////
376 // This is the geometry used for the ITS Pre-TDR and comes from an
377 // Euclid to Geant conversion. The only difference
378 // is in the details of the ITS supports. The detectors elements,
379 // detector numbering, and local and global reference frames are shown in
380 // the following figures.
383 <img src="picts/ITS/its1+2_convention_front_5.gif">
386 <font size=+2 color=red>
387 <p>This shows the front view of the SPDs.
390 <img src="picts/ITS/its1+2_convention_side_5.gif">
393 <font size=+2 color=red>
394 <p>This shows the perspective view of the SPDs.
396 <img src="picts/ITS/its1+2_tree.gif">
399 <font size=+2 color=red>
400 <p>This shows the geometry Tree for the SPDs.
405 <img src="picts/ITS/its3+4_convention_front_5.gif">
408 <font size=+2 color=red>
409 <p>This shows the front view of the SDDs.
412 <img src="picts/ITS/its3+4_convention_side_5.gif">
415 <font size=+2 color=red>
416 <p>This shows the perspective view of the SDDs.
418 <img src="picts/ITS/its3+4_tree.gif">
421 <font size=+2 color=red>
422 <p>This shows the geometry Tree for the SDDs.
427 <img src="picts/ITS/its5+6_convention_front_5.gif">
430 <font size=+2 color=red>
431 <p>This shows the front view of the SSDs.
434 <img src="picts/ITS/its5+6_convention_side_5.gif">
437 <font size=+2 color=red>
438 <p>This shows the perspective view of the SSDs.
441 <img src="picts/ITS/its5+6_tree.gif">
444 <font size=+2 color=red>
445 <p>This shows the geometry Tree for the SSDs.
450 <img src="picts/ITS/its_layer1-6_2.gif">
453 <font size=+2 color=red>
454 <p>This shows the front view of the whole ITS..
458 <img src="picts/ITS/its_layer1-6_1.gif">
461 <font size=+2 color=red>
462 <p>This shows the perspective view of the whole ITS..
466 <img src="picts/ITS/its1-6_tree.gif">
469 <font size=+2 color=red>
470 <p>This shows the geometry Tree for the whole ITS.
477 // Here are shown the details of the ITS support cones and services.
478 // First is a GEANT tree showing the organization of all of the volumes
479 // that make up the ITS supports and services.
482 <img src="picts/ITS/supports_tree.gif">
485 // What follows are a number of figures showing what these support
486 // structures look like.
490 <img src="picts/ITS/supports_3.gif">
493 <font size=+2 color=red>
494 <p>This shows the geometry of the supports for the Drift and Strip layers only.
498 <img src="picts/ITS/supports_2.gif">
501 <font size=+2 color=red>
502 <p>This shows the geometry of the supports for the Drift and Strip layers in front cut out.
506 <img src="picts/ITS/supports_1.gif">
509 <font size=+2 color=red>
510 <p>This shows the geometry of the supports for the Drift and Strip layers in a back cut out..
514 <img src="picts/ITS/suppssd.gif">
517 <font size=+2 color=red>
518 <p>This shows the geometry for the Strip layers supports.
522 <img src="picts/ITS/suppsdd.gif">
525 <font size=+2 color=red>
526 <p>This shows the geometry for the Drift layers supports.
532 // Read a file containing the geometry for the ITS version 5.
533 ////////////////////////////////////////////////////////////////////////
539 filtmp = gSystem->ExpandPathName(fEuclidGeometry.Data());
540 FILE *file = fopen(filtmp,"r");
544 printf("Ready to read Euclid geometry file\n");
545 ReadEuclid(fEuclidGeometry.Data(),topvol);
546 printf("Read in euclid geometries\n");
548 Error("CreateGeometry"," THE GEOM FILE %s DOES NOT EXIST !",
549 fEuclidGeometry.Data());
553 // Place the ITS ghost volume ITSV in its mother volume (ALIC) and make it
556 gMC->Gspos("ITSV",1,"ALIC",0,0,0,0,"ONLY");
558 // Outputs the geometry tree in the EUCLID/CAD format if requested to do so
561 gMC->WriteEuclid("ITSgeometry", "ITSV", 1, 5);
562 } // end if (fEuclidOut)
564 // read in the file containing the transformations for the active
565 // volumes for the ITS version 5. This is expected to be in a file
566 // ending in .det. This geometry is kept in the AliITSgeom class.
567 filtmp = gSystem->ExpandPathName(fEuclidGeometry.Data());
568 size = strlen(filtmp);
570 filtmp[size-3] = 'd'; // change from .euc to .det
571 filtmp[size-2] = 'e';
572 filtmp[size-1] = 't';
573 file = fopen(filtmp,"r");
574 if(file){ // if file exists use it to fill AliITSgeom structure.
576 printf("ready to read .det file %s\n",filtmp);
577 fITSgeom = new AliITSgeom(filtmp);
580 // fill AliITSgeom structure from geant structure just filled above
584 printf("finished with euclid geometrys\n");
586 //_____________________________________________________________________________
587 void AliITSv5symm::Init(){
588 ////////////////////////////////////////////////////////////////////////
589 // Initialise the ITS after it has been created.
590 ////////////////////////////////////////////////////////////////////////
596 //_____________________________________________________________________________
597 void AliITSv5symm::StepManager(){
598 ////////////////////////////////////////////////////////////////////////
599 // Called for every step in the ITS, then calles the AliITShit class
600 // creator with the information to be recoreded about that hit.
601 // The value of the macro ALIITSPRINTGEOM if set to 1 will allow the
602 // printing of information to a file which can be used to create a .det
603 // file read in by the routine CreateGeometry(). If set to 0 or any other
604 // value except 1, the default behavior, then no such file is created nor
605 // it the extra variables and the like used in the printing allocated.
606 ////////////////////////////////////////////////////////////////////////
611 TLorentzVector position, momentum;
612 TClonesArray &lhits = *fHits;
613 #if ALIITSPRINTGEOM==1
616 Float_t xl[3],xt[3],angl[6];
617 // Float_t par[20],att[20];
619 static Bool_t first=kTRUE,printit[6][50][50];
620 if(first){ for(copy1=0;copy1<6;copy1++)for(copy2=0;copy2<50;copy2++)
621 for(id=0;id<50;id++) printit[copy1][copy2][id] = kTRUE;
629 if(gMC->IsTrackInside()) vol[3] += 1;
630 if(gMC->IsTrackEntering()) vol[3] += 2;
631 if(gMC->IsTrackExiting()) vol[3] += 4;
632 if(gMC->IsTrackOut()) vol[3] += 8;
633 if(gMC->IsTrackDisappeared()) vol[3] += 16;
634 if(gMC->IsTrackStop()) vol[3] += 32;
635 if(gMC->IsTrackAlive()) vol[3] += 64;
637 // Fill hit structure.
638 if(!(gMC->TrackCharge())) return;
640 // Only entering charged tracks
641 if((id = gMC->CurrentVolID(copy)) == fIdSens[0]) {
643 id = gMC->CurrentVolOffID(0,copy);
644 //detector copy in the ladder = 1<->4 (ITS1)
646 gMC->CurrentVolOffID(1,copy1);
647 //ladder copy in the module = 1<->2 (I186)
648 gMC->CurrentVolOffID(2,copy2);
649 //module copy in the layer = 1<->10 (I132)
650 vol[2] = copy1+(copy2-1)*2;//# of ladders in one module = 2
651 } else if(id == fIdSens[1]){
653 id = gMC->CurrentVolOffID(0,copy);
654 //detector copy in the ladder = 1<->4 (ITS2)
656 gMC->CurrentVolOffID(1,copy1);
657 //ladder copy in the module = 1<->4 (I131)
658 gMC->CurrentVolOffID(2,copy2);
659 //module copy in the layer = 1<->10 (I132)
660 vol[2] = copy1+(copy2-1)*4;//# of ladders in one module = 4
661 } else if(id == fIdSens[2]){
663 id = gMC->CurrentVolOffID(1,copy);
664 //detector copy in the ladder = 1<->5 (ITS3 is inside I314)
666 id = gMC->CurrentVolOffID(2,copy);
667 //ladder copy in the layer = 1<->12 (I316)
669 } else if(id == fIdSens[3]){
671 id = gMC->CurrentVolOffID(1,copy);
672 //detector copy in the ladder = 1<->8 (ITS4 is inside I414)
674 id = gMC->CurrentVolOffID(2,copy);
675 //ladder copy in the layer = 1<->22 (I417)
677 }else if(id == fIdSens[4]){
679 id = gMC->CurrentVolOffID(1,copy);
680 //detector copy in the ladder = 1<->23 (ITS5 is inside I562)
682 id = gMC->CurrentVolOffID(2,copy);
683 //ladder copy in the layer = 1<->34 (I565)
685 }else if(id == fIdSens[5]){
687 id = gMC->CurrentVolOffID(1,copy);
688 //detector copy in the ladder = 1<->26 (ITS6 is inside I566)
690 id = gMC->CurrentVolOffID(2,copy);
691 //ladder copy in the layer = 1<->38 (I569)
694 return; // not an ITS volume?
695 } // end if/else if (gMC->CurentVolID(copy) == fIdSens[i])
697 gMC->TrackPosition(position);
698 gMC->TrackMomentum(momentum);
706 hits[7]=gMC->TrackTime();
707 // Fill hit structure with this new hit.
708 new(lhits[fNhits++]) AliITShit(fIshunt,gAlice->CurrentTrack(),vol,hits);
709 #if ALIITSPRINTGEOM==1
710 if(printit[vol[0]][vol[2]][vol[1]]){
711 printit[vol[0]][vol[2]][vol[1]] = kFALSE;
712 xl[0] = xl[1] = xl[2] = 0.0;
714 for(i=0;i<9;i++) mat[i] = 0.0;
715 mat[0] = mat[4] = mat[8] = 1.0; // default with identity matrix
718 gMC->Gdtom(xl,&(mat[0]),2);
721 gMC->Gdtom(xl,&(mat[3]),2);
724 gMC->Gdtom(xl,&(mat[6]),2);
726 angl[0] = TMath::ACos(mat[2]);
727 if(mat[2]==1.0) angl[0] = 0.0;
728 angl[1] = TMath::ATan2(mat[1],mat[0]);
729 if(angl[1]<0.0) angl[1] += 2.0*TMath::Pi();
731 angl[2] = TMath::ACos(mat[5]);
732 if(mat[5]==1.0) angl[2] = 0.0;
733 angl[3] = TMath::ATan2(mat[4],mat[3]);
734 if(angl[3]<0.0) angl[3] += 2.0*TMath::Pi();
736 angl[4] = TMath::ACos(mat[8]);
737 if(mat[8]==1.0) angl[4] = 0.0;
738 angl[5] = TMath::ATan2(mat[7],mat[6]);
739 if(angl[5]<0.0) angl[5] += 2.0*TMath::Pi();
741 for(i=0;i<6;i++) angl[i] *= 180.0/TMath::Pi(); // degrees
742 fp = fopen("ITSgeometry_v5.det","a");
743 fprintf(fp,"%2d %2d %2d %9e %9e %9e %9e %9e %9e %9e %9e %9e ",
744 vol[0],vol[2],vol[1], // layer ladder detector
745 xt[0],xt[1],xt[2], // Translation vector
746 angl[0],angl[1],angl[2],angl[3],angl[4],angl[5] // Geant rotaion
749 fprintf(fp,"%9e %9e %9e %9e %9e %9e %9e %9e %9e",
750 mat[0],mat[1],mat[2],mat[3],mat[4],mat[5],mat[6],mat[7],mat[8]
751 ); // Adding the rotation matrix.
754 } // end if printit[layer][ladder][detector]
759 //____________________________________________________________________________
760 void AliITSv5symm::Streamer(TBuffer &R__b){
761 ////////////////////////////////////////////////////////////////////////
762 // A dummy Streamer function for this class AliITSv5symm. By default it
763 // only streams the AliITS class as it is required. Since this class
764 // dosen't contain any "real" data to be saved, it doesn't.
765 ////////////////////////////////////////////////////////////////////////
767 if (R__b.IsReading()) {
768 Version_t R__v = R__b.ReadVersion();
770 AliITS::Streamer(R__b);
772 AliITS::Streamer(R__b);
775 R__b.WriteVersion(AliITSv5symm::IsA());
776 AliITS::Streamer(R__b);
777 } // end if R__b.IsReading()