]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSPlotGeom.C
Methos SetType added
[u/mrichter/AliRoot.git] / ITS / AliITSPlotGeom.C
1 void AliITSPlotGeom(const char *filename="galice.root"){
2 /////////////////////////////////////////////////////////////////////////
3 //   This macro displays part of the ITS sensitive volume as defined
4 // in the root file.
5 //   
6 //     Root > .L AliITSPlotGeom.C         //this loads the macro in memory
7 //     Root > AliITSPlotGeom();           //by default process first event   
8 //     Root > AliITSPlotGeom("galice2.root"); //process third event from 
9 //                                              galice2.root file.
10 //Begin_Html
11 /*
12 <img src="picts/AliITSplotgeom.gif">
13 */
14 //End_Html
15 /////////////////////////////////////////////////////////////////////////
16     if(gAlice){
17         delete gAlice;
18         gAlice=0;
19     }else{
20         // Dynamically link some shared libs
21         if(gClassTable->GetID("AliRun") < 0) {
22             gROOT->LoadMacro("loadlibs.C");
23             loadlibs();
24         } // end if
25     } // end if gAlice
26     // Connect the Root Galice file containing Geometry, Kine and Hits
27     TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(filename);
28     if(!file) file = new TFile(filename);
29
30     // Get AliRun object from file or create it if not on file
31     if(!gAlice) {
32         gAlice = (AliRun*)file->Get("gAlice");
33         if(gAlice) printf("AliRun object found on file\n");
34         if(!gAlice) gAlice = new AliRun("gAlice","Alice test program");
35     } // end if !gAlice
36       
37
38 // Get pointers to the ITS Alice detectors and Hits containers
39     AliITS    *ITS    = (AliITS*)    gAlice->GetDetector("ITS");
40     if(ITS==0){
41         cout << "ITS not found. Exiting." << endl;
42         return;
43     } // end if ITS==0
44     AliITSgeom *itsgeom = ITS->GetITSgeom();
45     if(itsgeom==0){
46         cout << "ITSgeom not defined. Exiting." << endl;
47         return;
48     } // end if itsgeom==0
49     Int_t is,ie,in;
50     is = itsgeom->GetStartSPD();
51     ie = itsgeom->GetLastSPD()+1;
52     in = ie-is;
53
54     // Transformations.
55 //    Float_t *x = new Float_t[in];
56 //    Float_t *y = new Float_t[in];
57 //    Float_t *z = new Float_t[in];
58     Float_t x,y,z;
59     TRotMatrix **r = new TRotMatrix*[in];
60     Int_t i,j;
61     Double_t m[9];
62     char name[10],title[20],name2[10],title2[20];
63 //    cout << in << endl;
64     for(i=0;i<in;i++){
65         itsgeom->GetRotMatrix(i+is,m);
66         sprintf(name,"ROT%d",i+is);
67         sprintf(title,"ROT%d",i+is);
68         r[i] = new TRotMatrix(name,title,m);
69 //      cout << name << title << endl;
70 //      itsgeom->GetTrans(i+is,x[i],y[i],z[i]);
71 //      cout << i << " " << x[i] << " " << y[i] << " " << z[i] <<endl;
72     } // end for i
73
74     // Sensitive volumes.
75     if(itsgeom->IsShapeDefined(0)){
76         TBRIK *spds = (TShape*) (itsgeom->GetShape(0));
77     }else{
78         TBRIK *spds = new TBRIK("SPD","SPD","void",0.64,0.015,3.48);
79     } // end if
80 //    cout << spds << endl;
81     if(itsgeom->IsShapeDefined(1)){
82         TBRIK *sdds = (TShape*) (itsgeom->GetShape(240));
83     }else{
84         TBRIK *sdds = new TBRIK("SDD","SDD","void",3.50850,0.01499,3.76320);
85     } // end if
86 //    cout << sdds << endl;
87     if(itsgeom->IsShapeDefined(2)){
88         TBRIK *ssds = (TShape*) (itsgeom->GetShape(1000));
89     }else{
90         TBRIK *ssds = new TBRIK("SSD","SSD","void",3.65,0.015,2.00);
91     } // end if
92 //    cout << ssds << endl;
93
94     // Set up display.
95
96     TCanvas *c1 = new TCanvas("c1","ITS geometry",10,10,700,700);
97 //  TPad *p1 = new TPad("p1","p1",0.01,0.01,0.99,0.99,46,3,1);
98 //  p1->Draw();
99 //  p1->cd();
100 //  TView *view = new TView(1);
101 //  view->SetRange(-5,-5,-5,5,5,5);
102       TShape  *mother = new TBRIK("Mother","Mother Volume","void",10,10,10);
103 //      TShape  *mother = new TBRIK("Mother","Mother Volume","void",30,30,30);
104 //      TShape  *mother = new TBRIK("Mother","Mother Volume","void",50,50,50);
105       TNode *node = new TNode("node","Mother Node",mother);
106       node->cd();
107
108     // Set up nodes
109     TNode **itsn = new TNode*[in];
110     TPolyLine3D **pl = new TPolyLine3D*[in];
111 /*    Double_t p[5*3],axis[5*3]={1.,0.,0.,
112                                0.,0.,0.,
113                                0.,1.,0.,
114                                0.,0.,0.,
115                                0.,0.,1.}; */
116     Double_t p[19*3],axis[19*3]={-0.25,0.,1.25,  // 1
117                                  +0.00,0.,1.25,  // 2
118                                  -0.25,0.,1.00,  // 3
119                                  +0.00,0.,1.00,  // 4
120                                  +0.00,0.,0.00,  // 5
121                                  +1.00,0.,0.00,  // 6
122                                  +1.25,0.,0.25,  // 7
123                                  +1.125,0.,0.125,  // 8
124                                  +1.00,0.,0.25,  // 9
125                                  +1.125,0.,0.125,  // 10
126                                  +1.25,0.,0.00,  // 11
127                                  +1.125,0.,0.125,  // 12
128                                  +1.00,0.,0.00,  // 13
129                                  +0.00,0.,0.00,  // 14
130                                  +0.00,1.,0.00,  // 15
131                                  +0.00,1.,0.125,  // 16
132                                  +0.25,1.,0.25,  // 17
133                                  +0.00,1.,0.125,  // 18
134                                  -0.25,1.,0.25,  // 19
135     };
136     for(i=0;i<in;i++){
137         itsgeom->GetTrans(i+is,x,y,z);
138 /*      switch (itsgeom->GetGeomMatrix(i+is)->GetDetectorIndex())
139         case 0:
140             sprintf(name,"SPD%d",i+is);
141             sprintf(title,"SPD%d",i+is);
142             sprintf(name2,"BSPD%d",i+is);
143             sprintf(title2,"BSPD%d",i+is);
144             itsn[i] = new TNode(name,title,new TBRIK(name2,title2,"void",
145                       spds->GetDx(),spds->GetDy(),spds->GetDz()),x,y,z,r[i]);
146             break;
147         case 1:
148             sprintf(name,"SDD%d",i+is);
149             sprintf(title,"SDD%d",i+is);
150             sprintf(name2,"BSDD%d",i+is);
151             sprintf(title2,"BSDD%d",i+is);
152             itsn[i] = new TNode(name,title,new TBRIK(name2,title2,"void",
153                       sdds->GetDx(),sdds->GetDy(),sdds->GetDz()),x,y,z,r[i]);
154             break;
155         case 2:
156             sprintf(name,"SSD%d",i+is);
157             sprintf(title,"SSD%d",i+is);
158             sprintf(name2,"BSSD%d",i+is);
159             sprintf(title2,"BSSD%d",i+is);
160             itsn[i] = new TNode(name,title,new TBRIK(name2,title2,"void",
161                       ssds->GetDx(),ssds->GetDy(),ssds->GetDz()),x,y,z,r[i]);
162             break; */
163         for(j=0;j<19;j++) itsgeom->LtoG(i+is,(Double_t*)&axis[3*j],
164                                             (Double_t*)&p[3*j]);
165         pl[i] = new TPolyLine3D(19,p);
166     } // end for i
167
168     // display it
169     node->cd();
170     node->Draw();
171 //    for(i=0;i<in;i++) itsn[i]->Draw();
172 //    node->Draw();
173     for(i=0;i<in;i++) pl[i]->Draw();
174     c1->Update();
175
176     // clean up
177 //    delete[] x;
178 //    delete[] y;
179 //    delete[] z;
180 //    for(i=0;i<itsgeom->GetIndexMax();i++) delete[] r[i];
181 //    delete[] r;
182 }
183