StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StFgtCosmicTrackPlots.cxx
1 /***************************************************************************
2  *
3  * $Id: StFgtCosmicTrackPlots.cxx,v 1.2 2012/01/31 23:37:17 avossen Exp $
4  * Author: C. K. Riley (ckriley@bnl.gov), Oct. 18 2011
5  *
6  ***************************************************************************
7  *
8  * Description: Plotting class for different histograms concerning
9  * cosmicstand tracks of middle quadrant
10  *
11  ***************************************************************************
12  *
13  * $Log: StFgtCosmicTrackPlots.cxx,v $
14  * Revision 1.2 2012/01/31 23:37:17 avossen
15  * fixed StFgtCosmicTrackMaker paths
16  *
17  * Revision 1.1 2012/01/31 23:25:35 avossen
18  * moved StFgtCosmicTrackMaker to StFgtPool
19  *
20  * Revision 1.3 2012/01/31 16:56:43 wwitzke
21  * Changing for cosmic test stand.
22  *
23  * Revision 1.2 2012/01/26 13:13:11 sgliske
24  * Updated to use StFgtConsts, which
25  * replaces StFgtEnums and StFgtGeomDefs
26  *
27  * Revision 1.1 2012/01/17 20:03:05 sgliske
28  * moved from StFgtQaMakers
29  *
30  * Revision 1.12 2011/12/07 17:19:59 ckriley
31  * minor update
32  *
33  * Revision 1.11 2011/11/25 20:20:04 ckriley
34  * added statusmaker functionality
35  *
36  * Revision 1.10 2011/11/09 21:03:19 ckriley
37  * working version with current containers
38  *
39  * Revision 1.9 2011/11/01 19:07:24 sgliske
40  * Updated to correspond with StEvent containers, take 2.
41  *
42  * Revision 1.7 2011/10/26 20:39:18 ckriley
43  * compilation fix
44  *
45  * Revision 1.6 2011/10/25 20:43:32 ckriley
46  * added more plots
47  *
48  * Revision 1.5 2011/10/25 15:38:51 ckriley
49  * more plots
50  *
51  * Revision 1.3 2011/10/20 17:13:44 ckriley
52  * major update -> headers, tracks stored in StFgtEvent instead of StFgtDisc, changes to trackmaker and algorithms
53  *
54  *
55  **************************************************************************/
56 
57 #include "StFgtCosmicTrackPlots.h"
58 
59 #include <string>
60 #include <TH1F.h>
61 #include <TProfile.h>
62 #include <TROOT.h>
63 #include <TStyle.h>
64 #include <TCanvas.h>
65 #include "StMaker.h"
66 
67 #include "StRoot/StFgtPool/StFgtCosmicTestStandGeom/StFgtCosmicTestStandGeom.h"
68 #include "StRoot/StFgtPool/StFgtCosmicTrackMaker/StFgtCosmicTrackMaker.h"
69 #include "StRoot/StFgtPool/StFgtCosmicTrackMaker/StFgtCosmicTrack.h"
70 
71 #include "StRoot/StEvent/StEvent.h"
72 #include "StRoot/StEvent/StFgtCollection.h"
73 #include "StRoot/StEvent/StFgtHit.h"
74 #include "StRoot/StEvent/StFgtPoint.h"
75 
76 // constructors
77 StFgtCosmicTrackPlots::StFgtCosmicTrackPlots( const Char_t* name,
78  const Char_t* cosmicTrackerName,
79  Short_t discId,
80  Short_t quadId,
81  const Char_t* quadName ) :
82  StMaker( name ), mFgtTrackMakerName( cosmicTrackerName ), mTrackVecPtr( 0 ),
83  mDiscId( discId ), mQuadId( quadId ),
84  mXbins( 200 ), mXmin( -50 ), mXmax( 50 ),
85  mBinFactorX( 1 ), mQuadName( quadName ),
86  mPath(""), mFileNameBase(""), mFileNameKey( quadName )/*,mHist(0)*/ {
87 
88  // set the style to plain
89  gROOT->SetStyle("Plain");
90 
91  // few defaults
92  mXbins = 60;
93  mXmin = -30;
94  mXmax = 30;
95 };
96 
97 StFgtCosmicTrackPlots::StFgtCosmicTrackPlots(const StFgtCosmicTrackPlots&){
98  std::cerr << "TODO" << endl;
99  throw 0;
100 };
101 
102 // deconstructor
103 StFgtCosmicTrackPlots::~StFgtCosmicTrackPlots(){
104 /*
105  if( mHistX ) {
106  for( Int_t i=0;i<4;i++) {
107  for( Int_t j=0;j<4;j++) {
108  delete mHistX[i][j];
109  }
110  }
111  //delete mHistX;
112  }
113  if( mHistY ) {
114  for( Int_t i=0;i<4;i++) {
115  for( Int_t j=0;j<4;j++) {
116  delete mHistY[i][j];
117  }
118  }
119  //delete mHistY;
120  }
121  for( Int_t i=0;i<2;i++) {
122  if( mHist1Dcluster[i] )
123  mHist1Dcluster[i]=0;
124  }
125  if( mHist1Dclusters )
126  delete mHist1Dclusters;
127  if( mHistTracks )
128  delete mHistTracks;
129  if( mHist2Dpoint )
130  delete mHist2Dpoint;
131  if( mProfileEfficiency )
132  delete mProfileEfficiency;
133  if( mHist2Dcluster )
134  delete mHist2Dcluster;
135  if( mProfileChi2 )
136  delete mProfileChi2;
137  if( mProfileX )
138  delete mProfileX;
139  if( mProfileY )
140  delete mProfileY;
141 */
142 };
143 
144 // equals operator
145 StFgtCosmicTrackPlots& StFgtCosmicTrackPlots::operator=(const StFgtCosmicTrackPlots&){
146  std::cerr << "TODO" << endl;
147  throw 0;
148 };
149 
150 // initialize CosmicTrackPlots
151 Int_t StFgtCosmicTrackPlots::Init(){
152  Int_t ierr = kStOk;
153 
154  // getting track data
155  TObject *dataMaker = GetMaker( mFgtTrackMakerName.data() );
156  if( !dataMaker ){
157  LOG_FATAL << "::Init() could not get pointer to a maker with name '" << mFgtTrackMakerName << "'" << endm;
158  ierr = kStFatal;
159  };
160 
161  mTrackVecPtr = 0;
162  StFgtCosmicTrackMaker *trackMaker = 0;
163  if( !ierr ){
164  trackMaker = static_cast< StFgtCosmicTrackMaker* >( dataMaker );
165  mTrackVecPtr = &(trackMaker->getCosmicTrackVec());
166  };
167 
168  if( !ierr ){
169 /*
170  // Take care of the histogram
171  if( mHistX ) {
172  for( Int_t i=0;i<4;i++) {
173  for( Int_t j=0;j<4;j++) {
174  delete[] mHistX[i][j];
175  }
176  }
177  //delete mHistX;
178  }
179  if( mHistY ) {
180  for( Int_t i=0;i<4;i++) {
181  for( Int_t j=0;j<4;j++) {
182  delete[] mHistY[i][j];
183  }
184  }
185  //delete mHistY;
186  }
187  for( Int_t i=0;i<2;i++) {
188  if( mHist1Dcluster[i] )
189  mHist1Dcluster[i]=0;
190  }
191  if( mHist1Dclusters )
192  delete mHist1Dclusters;
193  if( mHistTracks )
194  delete mHistTracks;
195  if( mHist2Dpoint )
196  delete mHist2Dpoint;
197  if( mProfileEfficiency )
198  delete mProfileEfficiency;
199  if( mHist2Dcluster )
200  delete mHist2Dcluster;
201  if( mProfileChi2 )
202  delete mProfileChi2;
203  if( mProfileX )
204  delete mProfileX;
205  if( mProfileY )
206  delete mProfileY;
207 */
208 
209  // naming the plots
210  std::string name, nameforhistX, nameforhistY, titleX, titleY;
211  {
212  std::stringstream ss;
213  ss << GetName() << "_" << mDiscId << "_" << mQuadId;
214  ss >> name;
215  };
216 
217  // easiest way to ensure stream is cleared is construct a new one
218  {
219  std::stringstream sx, sy, ss, snx, sny;
220  std::string profTitleX, profTitleY, nameforprofX, nameforprofY;
221  sx << "#deltax vs X for Quad " << mQuadName << ";X cm;#deltax cm";
222  sy << "#deltay vs Y for Quad " << mQuadName << ";Y cm;#deltay cm";
223  snx << "X" << name;
224  sny << "Y" << name;
225  snx >> nameforprofX;
226  sny >> nameforprofY;
227  profTitleX = sx.str();
228  profTitleY = sy.str();
229  mProfileX = new TProfile( nameforprofX.data(), profTitleX.data(), 8, 0., 40. );
230  mProfileY = new TProfile( nameforprofY.data(), profTitleY.data(), 8, 0., 40. );
231  }
232 
233  // naming 1D clusters by plane -> 50 channels per bin
234 {
235  std::string name1D, name1DR, name1DP, title1D, title1DR, title1DP;
236  {
237  std::stringstream ss, sr, sp;
238  ss << GetName() << "_" << mDiscId << "_" << mQuadId;
239  sr << ss.str() << "_R";
240  sp << ss.str() << "_P";
241  ss << "_C";
242  ss >> name1D;
243  sr >> name1DR;
244  sp >> name1DP;
245  };
246  // 1D cluster plot binning
247  mXmin = 0;
248  mXmax = 6000;
249  mXbins = 40;
250 
251  // easiest way to ensure stream is cleared is construct a new one
252  std::stringstream ss, sr ,sp;
253  ss << " Quad " << mQuadName << " clusters";
254  sr << ss.str() << " in R; adc sum";
255  sp << ss.str() << " in Phi; adc sum";
256  ss << " in R & Phi; adc sum";
257  title1D = ss.str();
258  title1DR = sr.str();
259  title1DP = sp.str();
260  //cout << "title = " << title1DR << endl;
261 
262  mHist1Dclusters = new TH1F( name1D.data(), title1D.data(), mXbins/mBinFactorX, mXmin, /*mXmax*/8000 );
263  mHist1Dcluster[0] = new TH1F( name1DR.data(), title1DR.data(), mXbins/mBinFactorX, mXmin, mXmax );
264  mHist1Dcluster[1] = new TH1F( name1DP.data(), title1DP.data(), mXbins/mBinFactorX, mXmin, mXmax );
265 }
266 
267  // naming mHistTracks, mHistEfficiency, mHist2Dcluster and mHistChi2
268  {
269  std::stringstream st, st2, sc, se, s2, snt, snt2, snc, sne, sn2, ss, ssn;
270  std::string nameT, nameT2, nameC, nameE, name2, nameS,
271  titleT, titleT2, titleC, titleE, title2, titleS;
272  st << "Reconstructed tracks on Quad " << mQuadName << "; cm; cm";
273  snt << name << "_tracks";
274  st2 << "Reconstructed registered tracks on Quad " << mQuadName << "; cm; cm";
275  snt2 << name << "_realTracks";
276  se << "Efficiency of Quad " << mQuadName << "; cm; cm";
277  sne << name << "_efficiency";
278  s2 << "Cluster points on Quad " << mQuadName
279  << " of registered tracks; cm; cm";
280  sn2 << name << "_realPoints";
281  sc << "Chi2 values on Quad " << mQuadName << "; cm; cm";
282 // sc << "Distance from track to hit on Quad " << mQuadName << "; cm; cm";
283  snc << name << "_chi2";
284  ss << "Cluster points on Quad " << mQuadName << "; cm; cm";
285  ssn << name << "_points";
286  titleT = st.str();
287  nameT = snt.str();
288  titleT2= st2.str();
289  nameT2 = snt2.str();
290  titleE = se.str();
291  nameE = sne.str();
292  title2 = s2.str();
293  name2 = sn2.str();
294  titleC = sc.str();
295  nameC = snc.str();
296  titleS = ss.str();
297  nameS = ssn.str();
298 
299  mXmin = 0;
300  mXmax = 40;
301  mXbins = 8;
302  mBinFactorX = 2;
303 
304  mHistTracks = new TH2F( nameT.data(), titleT.data(), mXbins, mXmin, mXmax, mXbins, mXmin, mXmax );
305  mHistRealTracks = new TH2F( nameT2.data(), titleT2.data(), mXbins,mXmin, mXmax, mXbins, mXmin, mXmax );
306  mHist2Dpoint = new TH2F( nameS.data(), titleS.data(), mXbins, mXmin, mXmax, mXbins, mXmin, mXmax );
307  mProfileEfficiency = new TProfile2D( nameE.data(), titleE.data(), mXbins/mBinFactorX, mXmin, mXmax, mXbins/mBinFactorX, mXmin, mXmax );
308  mHist2Dcluster = new TH2F( name2.data(), title2.data(), mXbins*5, mXmin, mXmax, mXbins*5, mXmin, mXmax );
309  mProfileChi2 = new TProfile2D( nameC.data(), titleC.data(), mXbins, mXmin, mXmax, mXbins, mXmin, mXmax );
310  }
311 
312  mXmin = -30;
313  mXmax = 30;
314  mXbins = 60;
315  mBinFactorX = 1;
316 
317  // naming APV histograms
318  for(int i=0;i<4;i++) {
319  for(int j=0;j<4;j++) {
320  std::stringstream sx, sy, ss, sn, snx, sny;
321  snx << "X" << name;
322  sny << "Y" << name;
323  sx << " #deltax:";
324  sy << " #deltay:";
325  ss << " Quad " << mQuadName << " ";
326  if (i==0 && j==0) {ss << "X:0-10 Y:0-10;"; sn << "_1";}
327  else if (i==0 && j==1) {ss << "X:0-10 Y:10-20"; sn << "_5";}
328  else if (i==0 && j==2) {ss << "X:0-10 Y:20-30"; sn << "_9";}
329  else if (i==0 && j==3) {ss << "X:0-10 Y:30-40"; sn << "_13";}
330  else if (i==1 && j==0) {ss << "X:10-20 Y:0-10"; sn << "_2";}
331  else if (i==1 && j==1) {ss << "X:10-20 Y:10-20"; sn << "_6";}
332  else if (i==1 && j==2) {ss << "X:10-20 Y:20-30"; sn << "-10";}
333  else if (i==1 && j==3) {ss << "X:10-20 Y:30-40"; sn << "_14";}
334  else if (i==2 && j==0) {ss << "X:20-30 Y:0-10"; sn << "_3";}
335  else if (i==2 && j==1) {ss << "X:20-30 Y:10-20"; sn << "_7";}
336  else if (i==2 && j==2) {ss << "X:20-30 Y:20-30"; sn << "_11";}
337  else if (i==2 && j==3) {ss << "X:20-30 Y:30-40"; sn << "_15";}
338  else if (i==3 && j==0) {ss << "X:30-40 Y:0-10"; sn << "_4";}
339  else if (i==3 && j==1) {ss << "X:30-40 Y:10-20"; sn << "_8";}
340  else if (i==3 && j==2) {ss << "X:30-40 Y:20-30"; sn << "_12";}
341  else if (i==3 && j==3) {ss << "X:30-40 Y:30-40"; sn << "_16";}
342 
343  snx << sn.str();
344  sny << sn.str();
345  snx >> nameforhistX;
346  sny >> nameforhistY;
347  sx << ss.str() << ";#deltax cm";
348  sy << ss.str() << ";#deltay cm";
349  titleX = sx.str();
350  titleY = sy.str();
351  //cout << "title = " << title << endl;
352 
353  mHistX[j][i] = new TH1F( nameforhistX.data(), titleX.data(), mXbins/mBinFactorX, mXmin, mXmax );
354  mHistY[j][i] = new TH1F( nameforhistY.data(), titleY.data(), mXbins/mBinFactorX, mXmin, mXmax );
355  }
356  }
357  };
358  return ierr;
359 };
360 
361 // fill histograms with track data
363  Int_t ierr = kStOk;
364 
365  StEvent* eventPtr = 0;
366  eventPtr = (StEvent*)GetInputDS("StEvent");
367 
368  if( !eventPtr ) {
369  LOG_ERROR << "Error getting pointer to StEvent from '" << ClassName() << "'" << endm;
370  ierr = kStErr;
371  };
372 
373  StFgtCollection *fgtCollectionPtr = 0;
374 
375  if( eventPtr ) {
376  fgtCollectionPtr=eventPtr->fgtCollection();
377  };
378 
379  if( !fgtCollectionPtr) {
380  LOG_ERROR << "Error getting pointer to StFgtCollection from '" << ClassName() << "'" << endm;
381  ierr = kStErr;
382  };
383 
384  if( !ierr ){
385  Bool_t isTrack = false;
386 
387  StFgtCosmicTrackVec::iterator trackIter = mTrackVecPtr->begin();
388  for( trackIter = mTrackVecPtr->begin(); trackIter != mTrackVecPtr->end(); ++trackIter ){
389 
390  // extract track info
391  Float_t a = trackIter->getLineParameterA();
392  Float_t b = trackIter->getLineParameterB();
393  Float_t x_0 = trackIter->getLineParameterX_0();
394  Float_t y_0 = trackIter->getLineParameterY_0();
395  Float_t dX = trackIter->getVarX();
396  Float_t dY = trackIter->getVarY();
397  Float_t chi2 = trackIter->getChiSquare();
398  Float_t hitX = trackIter->getHitX();
399  Float_t hitY = trackIter->getHitY();
400  isTrack = trackIter->getIsTrack();
401  Float_t Z;
402  if(mDiscId==0) Z = 30.48; // 12 inches in cm, top quadrant
403  else if(mDiscId==1) Z = 15.24; // 6 inches in cm, middle quadrant
404  else if(mDiscId==2) Z = 0.; // take origin as bottom quadrant
405  else return kStFatal; // no info for this disc
406  if(abs(a)>0.5 || abs(b)>0.5)
407  continue; // look at cosmictrack near perpendicular
408  Float_t x = a*Z+x_0; // where track line passes through quad
409  Float_t y = b*Z+y_0;
410 
411  // now bin accordingly
412  Int_t region=0;
413  if(x<40) {
414  if(y<40) region=16;
415  if(y<30) region=12;
416  if(y<20) region= 8;
417  if(y<10) region= 4;
418  }
419  if(x<30) {
420  if(y<40) region=15;
421  if(y<30) region=11;
422  if(y<20) region= 7;
423  if(y<10) region= 3;
424  }
425  if(x<20) {
426  if(y<40) region=14;
427  if(y<30) region=10;
428  if(y<20) region= 6;
429  if(y<10) region= 2;
430  }
431  if(x<10) {
432  if(y<40) region=13;
433  if(y<30) region= 9;
434  if(y<20) region= 5;
435  if(y<10) region= 1;
436  }
437  if(region==1) {mHistX[0][0]->Fill(dX); mHistY[0][0]->Fill(dY);}
438  else if(region==2) {mHistX[0][1]->Fill(dX); mHistY[0][1]->Fill(dY);}
439  else if(region==3) {mHistX[0][2]->Fill(dX); mHistY[0][2]->Fill(dY);}
440  else if(region==4) {mHistX[0][3]->Fill(dX); mHistY[0][3]->Fill(dY);}
441  else if(region==5) {mHistX[1][0]->Fill(dX); mHistY[1][0]->Fill(dY);}
442  else if(region==6) {mHistX[1][1]->Fill(dX); mHistY[1][1]->Fill(dY);}
443  else if(region==7) {mHistX[1][2]->Fill(dX); mHistY[1][2]->Fill(dY);}
444  else if(region==8) {mHistX[1][3]->Fill(dX); mHistY[1][3]->Fill(dY);}
445  else if(region==9) {mHistX[2][0]->Fill(dX); mHistY[2][0]->Fill(dY);}
446  else if(region==10) {mHistX[2][1]->Fill(dX); mHistY[2][1]->Fill(dY);}
447  else if(region==11) {mHistX[2][2]->Fill(dX); mHistY[2][2]->Fill(dY);}
448  else if(region==12) {mHistX[2][3]->Fill(dX); mHistY[2][3]->Fill(dY);}
449  else if(region==13) {mHistX[3][0]->Fill(dX); mHistY[3][0]->Fill(dY);}
450  else if(region==14) {mHistX[3][1]->Fill(dX); mHistY[3][1]->Fill(dY);}
451  else if(region==15) {mHistX[3][2]->Fill(dX); mHistY[3][2]->Fill(dY);}
452  else if(region==16) {mHistX[3][3]->Fill(dX); mHistY[3][3]->Fill(dY);}
453  else continue;
454 
455  // filling histograms
456  mProfileX->Fill(x,dX);
457  mProfileY->Fill(y,dY);
458 
459  if(mDiscId==1) { // for middle disc
460  mHistTracks->Fill(x,y);
461  mProfileChi2->Fill(x,y,chi2);
462 
463  mHist2Dpoint->Fill( hitX , hitY );
464  if(isTrack) {
465  mHist2Dcluster->Fill( hitX , hitY );
466  mHistRealTracks->Fill(x,y);
467  }
468  Int_t isTrackInt = 0;
469  if(isTrack) isTrackInt = 1;
470  mProfileEfficiency->Fill(x,y,isTrackInt);
471  }
472  }
473 
474  if( isTrack ) {
475  StFgtHitCollection *hitCollectionPtr = 0;
476  if( !ierr ){
477  hitCollectionPtr = fgtCollectionPtr->getHitCollection( mDiscId );
478  };
479 
480  if( hitCollectionPtr ){
481  const StSPtrVecFgtHit &hitVec = hitCollectionPtr->getHitVec();
482  StSPtrVecFgtHitConstIterator hitIter;
483 
484  Short_t adcSum=0;
485  //LOG_INFO << "We have " << hitVec.size() << " clusters " << endm;
486  for( hitIter = hitVec.begin(); hitIter != hitVec.end(); ++hitIter ){
487  Bool_t haveR=false;
488 
489  if((*hitIter)->charge() < 200 ) continue;
490  haveR=((*hitIter)->getLayer()=='R');
491 
492  if( haveR ) mHist1Dcluster[0]->Fill( (*hitIter)->charge() );
493  else mHist1Dcluster[1]->Fill( (*hitIter)->charge() );
494  adcSum += (*hitIter)->charge();
495  };
496  if(adcSum!=0) mHist1Dclusters->Fill( adcSum );
497  };
498  };
499  };
500 
501  return ierr;
502 };
503 
504 // make the histogram pretty and save it to a file, if given a filename base
506  //cout << "StFgtCosmicTrackPlots::Finish()" << endl;
507 
508 /* old style of efficiency (not per event)
509  for(int i=0;i<8;i++) {
510  for(int j=0;j<8;j++) {
511  Float_t con2Dpoint = mHist2Dpoint -> GetBinContent(i+1,j+1);
512  Float_t conTrack = mHistTracks -> GetBinContent(i+1,j+1);
513  Float_t x = i*5+1; // for correctly binning into efficiency profile
514  Float_t y = j*5+1; // mHist2Dpoint and tracks have 8x8 bins
515  if(con2Dpoint!=0) mProfileEfficiency->Fill(x,y,conTrack/con2Dpoint);
516  }
517  }
518 */
519  if( !mFileNameBase.empty() ){
520 
521  // filename
522  std::string filename;
523  if( !mPath.empty() )
524  filename = mPath + "/";
525  filename += mFileNameBase;
526  if( !mFileNameKey.empty() )
527  filename += std::string( "." ) + mFileNameKey;
528  filename += ".eps";
529 
530  // draw
531  gROOT->SetStyle("Plain");
532  TCanvas *canX = new TCanvas( "fgtTrackQAcanX", "fgtTrackQAcanX", 900, 400 );
533  gStyle->SetOptStat(0);
534  gStyle->SetOptDate(0);
535  Int_t pad=1;
536  canX->Divide(4,4);
537  for( Int_t i=3;i>=0;i--) {
538  for( Int_t j=0;j<4;j++) {
539  canX->cd(pad);
540  mHistX[i][j]->Draw();
541  pad++;
542  }
543  }
544 
545  TCanvas *canY = new TCanvas( "fgtTrackQAcanY", "fgtTrackQAcanY", 900, 400 );
546  gStyle->SetOptStat(0);
547  gStyle->SetOptDate(0);
548  pad=1;
549  canY->Divide(4,4);
550  for( Int_t i=3;i>=0;i--) {
551  for( Int_t j=0;j<4;j++) {
552  canY->cd(pad);
553  mHistY[i][j]->Draw();
554  pad++;
555  }
556  }
557 
558  TCanvas *can1Dcluster = new TCanvas( "fgtTrackQAcan1Dcluster", "fgtTrackQAcan1Dcluster", 900, 400 );
559  gStyle->SetOptStat(100000);
560  gStyle->SetOptDate(0);
561  can1Dcluster->Divide(2,1);
562  can1Dcluster->cd(1);
563  mHist1Dcluster[0]->Draw();
564  can1Dcluster->cd(2);
565  mHist1Dcluster[1]->Draw();
566 
567  TCanvas *can1Dclusters = new TCanvas( "fgtTrackQAcan1Dclusters", "fgtTrackQAcan1Dclusters", 900, 400 );
568  gStyle->SetOptStat(100000);
569  gStyle->SetOptDate(0);
570  mHist1Dclusters->Draw();
571 
572  TCanvas *canTrack = new TCanvas( "fgtTrackQATrack", "fgtTrackQATrack", 900, 400 );
573  gStyle->SetOptStat(0);
574  gStyle->SetOptDate(0);
575  gStyle->SetPalette(1);
576  mHistTracks->Draw("colz");
577 
578  TCanvas *canRealTrack = new TCanvas( "fgtTrackQARealTrack", "fgtTrackQARealTrack", 900, 400 );
579  gStyle->SetOptStat(0);
580  gStyle->SetOptDate(0);
581  gStyle->SetPalette(1);
582  mHistTracks->Draw("colz");
583 
584  TCanvas *canEfficiency = new TCanvas( "fgtTrackQAEf", "fgtTrackQAEf", 900, 400 );
585  gStyle->SetOptStat(0);
586  gStyle->SetOptDate(0);
587  gStyle->SetPalette(1);
588  mProfileEfficiency->Draw("colz");
589 
590  TCanvas *can2Dcluster = new TCanvas( "fgtTrackQApoint", "fgtTrackQApoint", 900, 400 );
591  gStyle->SetOptStat(0);
592  gStyle->SetOptDate(0);
593  gStyle->SetPalette(1);
594  mHist2Dcluster->Draw("colz");
595 
596  TCanvas *canChi = new TCanvas( "fgtTrackQAChi", "fgtTrackQAChi", 900, 400 );
597  gStyle->SetOptStat(0);
598  gStyle->SetOptDate(0);
599  gStyle->SetPalette(1);
600  mProfileChi2->Draw("colz");
601 
602  TCanvas *canX2 = new TCanvas( "fgtTrackQAcanX2", "fgtTrackQAcanX2", 900, 400 );
603  gStyle->SetOptStat(0);
604  gStyle->SetOptDate(0);
605  mProfileX->Draw();
606 
607  TCanvas *canY2 = new TCanvas( "fgtTrackQAcanY2", "fgtTrackQAcanY2", 900, 400 );
608  gStyle->SetOptStat(0);
609  gStyle->SetOptDate(0);
610  mProfileY->Draw();
611 
612 
613  canX->Print( filename.data() );
614  canY->Print( filename.data() );
615  can1Dcluster->Print( filename.data() );
616  can1Dclusters->Print( filename.data() );
617  canTrack->Print( filename.data() );
618  canRealTrack->Print( filename.data() );
619  canEfficiency->Print( filename.data() );
620  can2Dcluster->Print( filename.data() );
621  canChi->Print( filename.data() );
622  canX2->Print( filename.data() );
623  canY2->Print( filename.data() );
624  delete canX;
625  delete canY;
626  delete can1Dcluster;
627  delete can1Dclusters;
628  delete canTrack;
629  delete canRealTrack;
630  delete canEfficiency;
631  delete can2Dcluster;
632  delete canChi;
633  delete canX2;
634  delete canY2;
635  };
636 
637  return kStOk;
638 };
639 
640 void StFgtCosmicTrackPlots::Clear( Option_t *opt ) {
641 
642 };
643 
644 
645 ClassImp(StFgtCosmicTrackPlots);
void Clear(Option_t *opt="")
User defined functions.
virtual const char * GetName() const
special overload
Definition: StMaker.cxx:237
Definition: Stypes.h:44
Definition: Stypes.h:41