StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
HitsDraw.C
1 //*-- Author : Valery Fine(fine@bnl.gov) 05/05/99
2 //
3 // Test for the new St_Points3D class family
4 // Derived from HitsDraw.C macro
5 //
6 // begin_html
7 // Full text of this macro can be downloaded from the repository: <a href = "http://www.star.bnl.gov/cgi-bin/cvsweb.cgi/StRoot/macros/HitsDraw.C">HitsDraw.C</a>
8 // end_html
9 // $Id: HitsDraw.C,v 1.7 1999/11/11 20:28:24 fisyak Exp $
10 // $Log: HitsDraw.C,v $
11 // Revision 1.7 1999/11/11 20:28:24 fisyak
12 // Remove BFC.C
13 //
14 // Revision 1.6 1999/07/09 17:54:12 fine
15 // A example of selecting subdetectors and range was introduced
16 //
17 // Revision 1.5 1999/06/03 00:48:59 fine
18 // New version of HitsDraw
19 //
20 // Revision 1.4 1999/04/27 21:42:07 fine
21 // Web ref. have been added
22 //
23 // Revision 1.3 1999/04/27 21:32:03 fine
24 // some new comments and clean up
25 //
26 // Revision 1.2 1999/04/27 15:29:58 fine
27 // Make up plus Polyline was replaced by PolyMarker3D
28 //
29 // Revision 1.1 1999/04/23 22:42:59 fine
30 // How to draw Geometry and hits all together
31 //
32 // Forward declarations
33 class St_Node;
34 class St_NodeView;
35 class St_g2t_tpc_hit;
36 
37 // Global variables to access from an interactive session
38 St_Node *hall = 0;
39 St_NodeView *fullView = 0;
40 St_NodeView *shortView = 0;
41 St_NodeView *sensible = 0;
42 St_NodeView *singleSector = 0;
43 St_g2t_tpc_hit *points = 0;
44 //________________________________________________________________________________
45 // Subroutine to define the track color
46 Color_t trackColor(Long_t geantId)
47 {
48  Color_t color = kYellow;
49  switch (geantId) {
50  case 1:
51  case 2:
52  case 3:
53  color = kRed;
54  break;
55  case 4:
56  case 5:
57  case 6:
58  color = kGreen;
59  break;
60  case 7:
61  case 8:
62  case 9:
63  color = kYellow;
64  break;
65  default:
66  color = kBlue;
67  };
68  return color;
69 }
70 
71 //________________________________________________________________________________
72 // Subroutine to draw a separate "HELP" TCanvas object
73 void x3dHelpDraw(){
74  TCanvas *paveCanvas = new TCanvas("x3d","X3D Help",500,20,540,220);
75  TPaveText *label = new TPaveText(0,0,1,1);
76  const Char_t *myMessage[] = {"Plot 3D view"
77  ,"Note: Under Windows NT OpenGL will be used by default instead x3d for UNIX"
78  ,"Select x3d windows and type \"m\" in there to get x3d HELP windows"
79  ,"DON NOT FORGET to type \"q\" letter to terminate x3d widget"
80  ,"to continue this session"};
81 
82  Int_t lenMessage = 0;
83  lenMessage = sizeof(myMessage)/4;
84  Int_t j;
85  label->Draw();
86  for (j=lenMessage-1;j>0;j--) {
87  label->InsertText(myMessage[j]);
88  gPad->Update();
89  }
90  paveCanvas->Modified();
91  paveCanvas->Update();
92 }
93 
94 //________________________________________________________________________________
95 // Loading the share libraries
96 void Load() {
97  gSystem->Load("St_base");
98  gSystem->Load("xdf2root");
99  gSystem->Load("St_Tables");
100 }
101 //________________________________________________________________________________
102 // Drawing subroutine
103 void HitsDraw(){
104  TCanvas *m_TreeD = new TCanvas("STAF","Events",10,600,200,200);
105  Load();
106 // Draw STAR detector first
107 // TWebFile f("http://www.star.bnl.gov/~fine/star.root");
108 //_______________________________________
109 //
110 // reading STAR GEANT geometry database
111 //_______________________________________
112 
113  TFile f("/star/u2a/fine/WWW/star.root");
114  // read STAR geometry database remotely
115  TGeometry *star = f.Get("STAR");
116  if (!star) {
117  printf("Sorry, STAR was not found !\n");
118  return;
119  }
120 //_______________________________________
121 //
122 // Remove hall from the list of ROOT nodes
123 // to make it free of ROOT control
124 //_______________________________________
125  TList *listOfNode = gGeometry->GetListOfNodes();
126  St_Node *hall = listOfNode->First();
127  // Remove hall from the list of ROOT nodes to make it free of ROOT control
128  listOfNode->Remove(hall);
129  listOfNode->Remove(hall);
130 //_______________________________________
131 //
132 // Create an iterator to navigate STAR geometry
133 //_______________________________________
134  St_DataSetIter volume(hall);
135 
136 
137 //_______________________________________
138 //
139 // Read XDF file with some table information
140 //_______________________________________
141  // St_XDFFile file("/disk00000/star/auau200/hijing135/default/b0_3/year2a/hadronic_on/tss_dst/set184_01_48evts.p1.xdf");
142  St_XDFFile file("/disk00000/star/test/dev/trs_Linux/Fri/year_1b/set0352_01_35evts.dst.xdf");
143  // skip first record
144  St_DataSet *skip = file.NextEventGet();
145  if (skip) delete skip;
146  // Read the second "real" record
147  St_DataSet *muons = file.NextEventGet();
148  if (!muons) return;
149 //_______________________________________
150 //
151 // Find in there two kind of the tables we whan to use now
152 //_______________________________________
153  St_DataSetIter next(muons);
154  const Char_t *table = "event/geant/Event/g2t_tpc_hit";
155  const Char_t *trackTable = "event/geant/Event/g2t_track";
156 
157 // struct g2t_tpc_hit_st {
158 // long id; // primary key
159 // long next_tr_hit_p; // Id of next hit on same track
160 // long track_p; // Id of parent track
161 // long volume_id; // STAR volume identification
162 // float de; // energy deposition at hit
163 // float ds; // path length within padrow
164 // float p[3]; // local momentum
165 // float tof; // time of flight
166 // float x[3]; // coordinate (Cartesian)
167 // }
168  points = (St_g2t_tpc_hit *)next(table);
169  if (!points) {
170  next.Du();
171  return;
172  }
173  points->Print(0,10);
174  g2t_tpc_hit_st *p = 0;
175  if (points)
176  g2t_tpc_hit_st *p = points->GetTable();
177  else {
178  printf("Error. The table <%s> was not found\n",table);
179  return;
180  }
181 
182 // struct g2t_track_st {
183 // . . . .
184 // long ge_pid; /* GEANT particle id */
185 // . . . .
186 //} G2T_TRACK_ST;
187 
188  St_g2t_track *trackT = (St_g2t_track *)next(trackTable);
189  trackT->Print(0,10);
190  g2t_track_st *geTrack = trackT->GetTable();
191 
192 //_______________________________________
193 //
194 // The table was found, let's create 3D viewer for that
195 //
196 // We want to see all x[0],x[1],x[2] columns of each row
197 // with one and the same value of "track_id" columns as one
198 // 3D objects.
199 //
200 // One objects means all part (points) of this simgle object
201 // will have one and same 3D attributes: color, size, style etc
202 //_______________________________________
203  St_Table3Points *track = 0;
204  TString tr;
205  tr = "track_p";
206  St_Table &ttt = *((St_Table *)points);
207  // Track2Line MUST be on heap othwesie 3D view will crash just code leaves this
208  // subroutine
209  St_TableSorter *Track2Line = new St_TableSorter (ttt,tr);
210 
211  Int_t i = 0;
212  Char_t buffer[10];
213  Int_t ntracks = 0;
214 // const Int_t maxtracks = 30;
215  const Int_t maxtracks = 5;
216 //---------------------------- Fill tracks -------------------
217  gBenchmark->Start("Fill time");
218  // Get Cave
219  long currentId = -1;
220  long newId = 0;
221  g2t_tpc_hit_st *hitPoint = 0;
222  printf(" new track found:");
223  St_Node *thisTrack[7] = {0,0,0,0,0,0,0}; // seven noded for 7 different colors
224  Int_t MaxRowtoCount = 5000; // 5000;
225  Int_t MaxTracks = Track2Line->CountKeys();
226  printf(" Total tracks here %d \n",MaxTracks);
227  MaxTracks = 100;
228 // Create one "dummy" node to hold all tracks
229  St_Node *allTracks = new St_Node(".hits","Global Hits", (TShape *)0);
230  if (!allTracks) { printf("Bug !!!\n"); return; }
231  allTracks->Mark();
232  allTracks->SetVisibility();
233  hall->Add(allTracks);
234  for (i=0;i<Track2Line->GetNRows() && ntracks < MaxTracks ;i++)
235  {
236  hitPoint = p + Track2Line->GetIndex(i);
237  newId = hitPoint->track_p;
238 
239  if (newId != currentId) {
240  printf(".");
241  const Char_t *xName = "x[0]";
242  const Char_t *yName = "x[1]";
243  const Char_t *zName = "x[2]";
244  track = new St_Table3Points(Track2Line,(const void *)&newId,xName,yName,zName);
245 
246  // Create a shape for this node
247  St_PolyLineShape *trackShape = new St_PolyLineShape(track);
248  trackShape->SetVisibility(1);
249  Int_t colorIndx = ntracks%7;
250 // trackShape->SetLineColor(trackColor(geTrack[newId]->ge_pid));
251  trackShape->SetColorAttribute(colorIndx+kGreen);
252  trackShape->SetLineStyle(1); trackShape->SetSizeAttribute(2);
253  // Create a node to hold it
254  if (!thisTrack[colorIndx]) {
255  thisTrack[colorIndx] = new St_Node("hits","hits",trackShape);
256  thisTrack[colorIndx]->Mark(); thisTrack[colorIndx]->SetVisibility();
257  St_NodePosition *pp = allTracks->Add(thisTrack[colorIndx]);
258  if (!pp) printf(" no position %d\n",ntrack);
259  }
260  else
261  thisTrack[colorIndx]->Add(trackShape);
262  currentId = newId;
263  ntracks++;
264  }
265  }
266  printf("\n");
267  printf("%d tracks have been read\n",ntracks);
268  gBenchmark->Stop("Fill time");
269 
270  printf(" Now we will try to draw something \n");
271  Int_t irep;
272  gBenchmark->Stop("Draw time");
273 // Select sensors
274  ((St_Node *)volume.FindByName("ZCAL"))->SetVisibility(0);
275  ((St_Node *)volume.FindByName("QDIV"))->SetVisibility(0);
276  fullView = new St_NodeView(*hall);
277  // Create the "open" sub-structure from the full one
278  sensible = new St_NodeView(fullView);
279  printf(" drawing the STAR geometry sensible volumes and hits \n");
280  sensible->Draw();
281 
282  // Begin_Html <P ALIGN=CENTER> <IMG SRC="gif/HitsDrawFullView.gif"> </P> End_Html //
283 
284  gPad->Update();
285  // Draw a help for x3d picture
286  printf(" drawing the HELP windows\n");
287 // x3dHelpDraw();
288  // make simple view
289 
290  // Select node to be left
291  printf(" Marking the current structure to create another simplified one\n");
292  sensible->Mark(); // Select HALL
293  shortView = (St_NodeView *)sensible->FindByName("TPGV"); // Select TPGV
294  if (shortView) shortView->Mark();
295 
296  // Select all hits // Select ALL Hits
297  sensible->FindByName(".hits")->Mark();
298  St_DataSetIter nextHits(sensible->FindByName(".hits"));
299  while (shortView = (St_NodeView *) nextHits()) shortView->Mark();
300 
301  nextHits.Reset(sensible);
302  while (shortView = (St_NodeView *) nextHits()) {
303  if (strcmp(shortView->GetName(),"ZCAL")==0) continue; // skip ZCAL detector element
304  shortView->Mark();
305  }
306 
307  // Create new short view // New "short" dataset
308  printf(" Creating a new structure simplified structure\n");
309  shortView = new St_NodeView(sensible);
310  TCanvas *m_simpleD = new TCanvas("Detector","Simple view and hits",500,500,400,400);
311  shortView->Draw("L");
312 
313  // Begin_Html <P ALIGN=CENTER> <IMG SRC="gif/HitsDrawSimpleView.gif"> </P> End_Html //
314 
315  printf(" Drawing the new structure simplified structure\n");
316  gPad->Update();
317  gBenchmark->Summary();
318  printf(" Drawing 3d solid view of the last structure\n");
319 //== m_simpleD->x3d();
320 
321  // Begin_Html <P ALIGN=CENTER> <IMG SRC="gif/HitsDrawX3D.gif"> </P> End_Html //
322 
323  printf(" The following global variables are available:\n");
324  printf("\tSt_Node *hall - \"closed\" full GEANT geometry structure\n");
325  printf("\tSt_NodeView *fullView - \"open\" full GEANT geometry structure\n");
326  printf("\tSt_NodeView *sensible - all \"sensible\" detector elements of the full GEANT structure\n");
327  printf("\tSt_NodeView *shortView - subset of the structure above\n");
328 
329  // Draw hits and one sector only
330  St_NodeView *singleSector = new St_NodeView(fullView,".hits","TPSS");
331  St_NodeView *sector = (St_NodeView *)singleSector->FindByName("TPSS");
332  sector->SetVisibility(2); // This node is invisibible but sons
333  TCanvas *ccc = new TCanvas("STARTPC","Sectors",500,10,200,200);
334  ccc.cd();
335  singleSector->Draw();
336  gPad->Update();
337  // Set a new range
338 #if 1
339  Float_t min[3];
340  Float_t max[3];
341  TVirtualPad *thisPad = gPad;
342  TView *view = thisPad->GetView();
343  sector->GetGlobalRange(singleSector,min,max);
344  view->SetRange(min,max);
345  // Update the last picture within the range provided
346  thisPad->Modified();
347  thisPad->Update();
348 #endif
349 
350 }
virtual TDataSet * First() const
Return the first object in the list. Returns 0 when list is empty.
Definition: TDataSet.cxx:403
virtual Size_t SetSizeAttribute(Size_t size)
to be documented
Int_t GetIndex(UInt_t sortedIndex) const
returns the original index of the row by its sorted index
virtual void Remove(TDataSet *set)
Remiove the &quot;set&quot; from this TDataSet.
Definition: TDataSet.cxx:641
virtual void SetVisibility(ENodeSEEN vis=TVolume::kBothVisible)
Definition: TVolume.cxx:753
virtual TDataSet * FindByName(const Char_t *name, const Char_t *path="", Option_t *opt="")
to be documented
virtual Color_t SetColorAttribute(Color_t color)
to be documented
Definition: TTable.h:48
virtual TDataSet * FindByName(const char *name, const char *path="", Option_t *opt="") const
Definition: TDataSet.cxx:378
virtual void Draw(Option_t *depth="3")
Draw Referenced node with current parameters.
virtual Int_t GetGlobalRange(const TVolumeView *rootNode, Float_t *min, Float_t *max)
virtual void SetVisibility(Int_t vis=1)
to be documented
virtual Int_t CountKeys() const