StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
tpcdraw.C
1 // $Id: tpcdraw.C,v 1.14 2000/03/08 22:52:49 snelling Exp $
2 // $Log: tpcdraw.C,v $
3 // Revision 1.14 2000/03/08 22:52:49 snelling
4 // Changed View settings (thanks Ian)
5 //
6 // Revision 1.13 2000/01/20 02:18:28 snelling
7 // fixed St_TableSorter under Linux still crash under SUN (CINT related)
8 //=======================================================================
9 // owner: Raimond Snellings
10 // what it does: Macro for plotting hits and pixels in combination with bfc.C
11 //=======================================================================
12 
13 //int draw_sector(TPad &padname, Float_t theta, Float_t phi) {
14 int InitDraw() {
15  /*
16  * -----------------------------------------------------------------
17  * draw TPC sectors, two sections, inner and outer.
18  * -----------------------------------------------------------------*/
19  Float_t RminInner = 51.9/2;
20  Float_t RmaxInner = 126.5/2;
21  Float_t RminOuter = 127.5/2;
22  Float_t RmaxOuter = 202.9/2;
23  Float_t Zmin = -210.;
24  Float_t Zmax = 210.;
25 
26  // Define the TPC sector outlines.
27  TPGON *InnerSector = new TPGON("InnerSector","TPC inner sector","void",-15,360,12,2);
28  InnerSector->DefineSection(0, Zmin, RminInner, RmaxInner);
29  InnerSector->DefineSection(1, Zmax, RminInner, RmaxInner);
30  InnerSector->SetLineColor(2);
31  TNode *m_node1 = new TNode("m_node1","inner sector",InnerSector);
32 
33  TPGON *OuterSector = new TPGON("OuterSector","TPC outer sector","void",-15,360,12,2);
34  OuterSector->DefineSection(0,Zmin,RminOuter,RmaxOuter);
35  OuterSector->DefineSection(1,Zmax,RminOuter,RmaxOuter);
36  OuterSector->SetLineColor(4);
37  TNode *m_node2 = new TNode("m_node2","outer sector",OuterSector);
38 
39  return kStOK;
40 }
41 //______________________________________________________________________
42 
43 int DrawEvent(TPad *padname, Float_t theta, Float_t phi) {
44 
45  // * Plot hits and tracks in the TPC */
46 
47  char tkstring[10];
48  Float_t zmin=-210;
49  Float_t zmax=210;
50  Float_t rmin[3]={ -200., -200., -200.};
51  Float_t rmax[3]={ 200., 200., 200.};
52 
53  TView *view = new TView(rmin,rmax);
54  padname->SetView(view);
55  view->RotateView(phi,theta);
56  // view->ShowAxis();
57  view->Draw();
58 
59  m_node2->Draw("same");
60  m_node1->Draw("same");
61 
62  // Create iterators for the two datasets
63  St_DataSetIter tpc_data(chain->DataSet("tpc_hits"));
64  St_tcl_tphit *hits = NULL;
65  tcl_tphit_st *hit1 = NULL;
66  hits = (St_tcl_tphit* ) tpc_data.Find("tphit");
67  if (hits) {hit1 = hits->GetTable();}
68  else { cout << "Error: tphit table header does not exist " << endl; return kStWarn; }
69  if (!hit1) { cout << "Error: tphit table does not exist " << endl; return kStWarn; }
70  Int_t nhits = hits->GetNRows();
71  if (nhits == 0) {cout << "Error: tphit table contains zero rows " << endl; return kStWarn;}
72  St_TableSorter sortrk(hits,"track",0,nhits-1);
73 
74  St_DataSetIter tpc_tracks(chain->DataSet("tpc_tracks"));
75  St_tpt_track *track = 0;
76  tpt_track_st *track1 = 0;
77  track = (St_tpt_track *) tpc_tracks.Find("tptrack");
78  if (track) {track1 = track->GetTable();}
79  else { cout << "Error: tptrack table header does not exist " << endl; return kStWarn; }
80  if (!track1) { cout << "Error: tptrack table does not exist " << endl; return kStWarn; }
81  Int_t ntracks = track->GetNRows();
82  if (ntracks == 0) {cout << "Error: tptrack table contains zero rows " << endl; return kStWarn;}
83 
84 
85  // Draw hits
86  Int_t k=0;
87  Float_t *pts = new Float_t[3*nhits];
88  for (Int_t j = 0; j < nhits; j++) {
89  if (hit1[j]->z > zmin && hit1[j]->z < zmax) {
90  pts[3*k] = hit1[j]->x;
91  pts[3*k+1] = hit1[j]->y;
92  pts[3*k+2] = hit1[j]->z;
93  k++;
94  }
95  }
96  cout << "total hits " << k << endl;
97  TPolyMarker3D *hit = new TPolyMarker3D(k,pts,2);
98  hit->SetMarkerColor(1);
99  hit->Draw("same");
100 
101  //Draw Hits on track with different colors
102  Int_t ngoodtrks = 0;
103  for (j = 0; j < ntracks; j++) {
104  Int_t trkid = track1[j]->id;
105  Int_t trkflag = track1[j]->flag;
106  Int_t trknrec = track1[j]->nrec;
107  cout <<" Track #" << j << " "<< trknrec <<" hits" << endl;
108  Float_t *tkpts = new Float_t[3*trknrec];
109  // define colors and markers for different tracks
110  Int_t mark = j + 20;
111  if (j>11) {mark = j + 8;}
112  if (j>22) {mark = j - 3;}
113  if (j>33) {mark = j - 14;}
114  if (j>44) {mark = j - 25;}
115  Int_t incolor = j + 2;
116  // colors 0 and 10 are white - invisible
117  if (incolor == 10) {incolor = 2;}
118  if (incolor > 50) {incolor -= 50;}
119  // loop over hits belonging to a track
120  int k =0;
121  for (Int_t i = 0; i < trknrec; i++) {
122  Int_t hitid = 1000*trkid + i + 1; // because table start at 1
123  Int_t irow_hit = sortrk[hitid];
124  if (irow_hit < 0) continue;
125  if (hit1[irow_hit]->z > zmin && hit1[irow_hit]->z < zmax) {
126  tkpts[3*k] = hit1[irow_hit]->x;
127  tkpts[3*k+1] = hit1[irow_hit]->y;
128  tkpts[3*k+2] = hit1[irow_hit]->z;
129  k++;
130  }
131  }
132  TPolyMarker3D *thit = new TPolyMarker3D(k,tkpts,mark);
133  thit->SetMarkerColor(incolor);
134  thit->Draw("same");
135 
136 
137  // Draw good tracks as polylines.
138  Float_t p[6];
139  if (trkflag < 0) {continue;}
140  Float_t ay = 1.0/tan(0.0174533*track1[j]->psi);
141  Float_t x0 = track1[j]->r0*cos(0.0174533*track1[j]->phi0);
142  Float_t y0 = track1[j]->r0*sin(0.0174533*track1[j]->phi0);
143  if (track1[j]->z0 >zmin && track1[j]->z0 < zmax) {
144  p[1] = rmin[1];
145  p[4] = rmax[1];
146  p[0] = x0 + ay*(p[1]-y0);
147  if (p[0] < rmin[0]) {
148  p[0] = rmin[0];
149  p[1] = y0 + (p[0]-x0)/ay;
150  }
151  if (p[0] > rmax[0]) {
152  p[0] = rmax[0];
153  p[1] = y0 + (p[0]-x0)/ay;
154  }
155  p[2] = track1[j]->tanl*(y0-p[1]) + track1[j]->z0;
156  p[3] = x0 + ay*(p[4]-y0);
157  if (p[3] < rmin[0]) {
158  p[3] = rmin[0];
159  p[4] = y0 + (p[3]-x0)/ay;
160  }
161  if (p[3] > rmax[0]) {
162  p[3] = rmax[0];
163  p[4] = y0 + (p[3]-x0)/ay;
164  }
165  p[5] = track1[j]->tanl*(y0-p[4]) + track1[j]->z0;
166  TPolyLine3D *trk = new TPolyLine3D(2,p);
167  trk->SetLineColor(incolor);
168  trk->Draw("same");
169  ngoodtrks++;
170  }
171  }
172 
173 
174  padname->Modified();
175  padname->Update();
176 
177  return kStOK;
178 }
179 
180 //-----------------------------------------------------------------------
181 
182 int DrawPixels(Text_t* varexp, Text_t* selection, Text_t* options) {
183  /*
184  * -----------------------------------------------------------------
185  * draw TPC pixels
186  * -----------------------------------------------------------------*/
187 
188  St_DataSetIter Itpc_raw(chain->DataSet("tpc_raw"));
189  // St_DataSetIter Itpc_raw(GetDataSet("tpc_raw"));
190  St_tfc_adcxyz *phtfc = 0;
191  tfc_adcxyz_st *ptadcxyz = 0;
192  phtfc = (St_tfc_adcxyz *) Itpc_raw.Find("adcxyz");
193 
194  if (phtfc) {ptadcxyz = phtfc->GetTable();}
195  else { cout << "Warning: adcxyz table header does not exist " << endl; return kStWarn; }
196  if (!ptadcxyz) { cout << "Warning: adcxyz table does not exist " << endl; return kStWarn; }
197 
198  TH1* pPixelHist = phtfc->Draw(varexp,selection,options);
199  pPixelHist->SetMarkerStyle(26);
200  pPixelHist->SetMarkerColor(4);
201 
202  return kStOK;
203 }
204 
205 //-----------------------------------------------------------------------
206 
207 int DrawHits(Text_t* varexp, Text_t* selection, Text_t* options){
208  /*
209  * -----------------------------------------------------------------
210  * draw TPC hits
211  * -----------------------------------------------------------------*/
212 
213  // get pointers to tpc hit table
214  St_DataSetIter Itpc_hits(chain->DataSet("tpc_hits"));
215  St_tcl_tphit *phtcl = 0;
216  tcl_tphit_st *pttphit = 0;
217  phtcl = (St_tcl_tphit *) Itpc_hits.Find("tphit");
218  if (phtcl) {pttphit = phtcl->GetTable();}
219  else { cout << "Error: tphit table header does not exist " << endl; return kStWarn; }
220  if (!pttphit) { cout << "Error: tphit table does not exist " << endl; return kStWarn; }
221 
222  TH1* pHitHist = phtcl->Draw(varexp,selection,options);
223  pHitHist->SetMarkerStyle(20);
224  pHitHist->SetMarkerColor(2);
225 
226  return kStOK;
227 }
228 
229 //-----------------------------------------------------------------------
230 
231 int DrawGeantHits(Text_t* varexp, Text_t* selection, Text_t* options) {
232  /*
233  * -----------------------------------------------------------------
234  * draw TPC geant hits
235  * -----------------------------------------------------------------*/
236 
237  St_DataSetIter Igeant_g2t(chain->DataSet("g2t_tpc_hit"));
238  St_g2t_tpc_hit *phg2t = 0;
239  g2t_tpc_hit_st *ptg2t_tpc_hit = 0;
240  phg2t = (St_g2t_tpc_hit *) Igeant_g2t.Find("g2t_tpc_hit");
241  if (phg2t) {ptg2t_tpc_hit = phg2t->GetTable();}
242  else { cout << "Warning: g2t tpc hit table header does not exist " << endl; return kStWarn; }
243  if (!ptg2t_tpc_hit) { cout << "Warning: g2t tpc hit table does not exist " << endl; return kStWarn; }
244 
245  TH1* pGeantHist = phg2t->Draw(varexp,selection,options);
246  pGeantHist->SetMarkerStyle(24);
247  pGeantHist->SetMarkerColor(3);
248 
249  return kStOK;
250 }
251 
252 //-----------------------------------------------------------------------
253 
254 void tpcdraw() {
255  /*
256  * -----------------------------------------------------------------
257  * draw TPC sectors with some options
258  * -----------------------------------------------------------------*/
259 
260  /*
261  if (chain->GetOption("miniDAQ")) {
262  chain->SetInput("BEGIN_RUN",".make/xdfin/.const/BEGIN_RUN");
263  chain->SetInput("ChainFlags[kTPC]_DATA",".make/xdfin/.data/ChainFlags[kTPC]_DATA");
264  }
265  */
266 
267  gStyle->SetCanvasColor(10); // white
268  gStyle->SetPadColor(10); // white
269  gStyle->SetFrameFillColor(10); // white
270  gStyle->SetOptFit(1);
271 
272  TCanvas *c1 = new TCanvas("tpcDraw","x y z info",100,100,700,700);
273  c1->SetFillColor(0);
274  c1->SetBorderMode(0);
275  c1->SetBorderSize(0);
276 
277  TPaveText *title = new TPaveText(.2,0.96,.8,.995);
278  title->SetFillColor(33);
279  title->AddText("pixels and hits");
280  title->Draw();
281 
282  TPad *pad1 = new TPad("pad1","x vs y",0.01,0.48,0.49,0.95);
283  TPad *pad2 = new TPad("pad2","x vs z",0.51,0.48,0.98,0.95);
284  TPad *pad3 = new TPad("pad3","tpc view 1",0.01,0.01,0.49,0.48);
285  TPad *pad4 = new TPad("pad4","tpc view 2",0.51,0.01,0.98,0.48);
286 
287  pad1->SetFillColor(0);
288  pad1->SetBorderMode(0);
289  pad1->SetBorderSize(0);
290  pad1->Draw();
291 
292  pad2->SetFillColor(0);
293  pad2->SetBorderMode(0);
294  pad2->SetBorderSize(0);
295  pad2->Draw();
296 
297  pad3->SetFillColor(0);
298  pad3->SetBorderMode(0);
299  pad3->SetBorderSize(0);
300  pad3->Draw();
301 
302  pad4->SetFillColor(0);
303  pad4->SetBorderMode(0);
304  pad4->SetBorderSize(0);
305  pad4->Draw();
306 
307  pad1->cd();
308  if (chain->GetOption("trs") ||
309  chain->GetOption("miniDAQ") || chain->GetOption("TDAQ")) {
310  DrawPixels("x:y","(adc>2)*adc","box");
311  }
312  if (DrawPixels == 0) { DrawHits("y:x","","same,scat"); }
313  else { DrawHits("x:y","","scat"); }
314  if (!chain->GetOption("miniDAQ") && !chain->GetOption("TDAQ")) {
315  DrawGeantHits("x[0]:x[1]","volume_id<2446","same,scat");
316  }
317 
318  pad2->cd();
319  if (chain->GetOption("trs") ||
320  chain->GetOption("miniDAQ") || chain->GetOption("TDAQ")) {
321  DrawPixels("x:z","(adc>2)*adc","box");
322  }
323  if (DrawPixels == 0) { DrawHits("z:x","","same,scat"); }
324  else { DrawHits("x:z","","scat"); }
325  if (!chain->GetOption("miniDAQ") && !chain->GetOption("TDAQ")) {
326  DrawGeantHits("x[0]:x[2]","volume_id<2446","same,scat");
327  }
328 
329  InitDraw();
330 
331  pad3->cd();
332  DrawEvent(pad3,0.,-90.);
333 
334  pad4->cd();
335  DrawEvent(pad4,90.,-90.);
336 }
337 
Definition: Stypes.h:42
Definition: Stypes.h:40
virtual TDataSet * Find(const Char_t *path, TDataSet *rootset=0, Bool_t mkdir=kFALSE, Bool_t titleFlag=kFALSE)