StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StHistUtil.cxx
1 // $Id: StHistUtil.cxx,v 2.110 2021/03/19 19:27:35 genevb Exp $
2 // $Log: StHistUtil.cxx,v $
3 // Revision 2.110 2021/03/19 19:27:35 genevb
4 // Small adjustment to label size
5 //
6 // Revision 2.109 2021/03/18 21:37:39 genevb
7 // Update styles of numerous histograms
8 //
9 // Revision 2.108 2019/05/22 21:24:30 genevb
10 // Add sDCA vs. time-in-run
11 //
12 // Revision 2.107 2019/03/26 15:29:35 genevb
13 // Introduce ETOF
14 //
15 // Revision 2.106 2019/03/14 02:31:52 genevb
16 // Introduce iTPC plots
17 //
18 // Revision 2.105 2019/03/04 20:55:59 genevb
19 // Improve RDO layout for iTPC
20 //
21 // Revision 2.104 2019/02/25 19:20:18 genevb
22 // Sector numbering fix
23 //
24 // Revision 2.103 2018/07/06 22:13:04 smirnovd
25 // [Cosmetic] Changes in white space
26 //
27 // Revision 2.102 2018/07/06 22:10:26 smirnovd
28 // [Cosmetic] Inverse test conditions to skip loop iterations
29 //
30 // Revision 2.101 2018/05/05 04:00:54 genevb
31 // iTPC RDO outlines
32 //
33 // Revision 2.100 2018/05/02 21:06:32 genevb
34 // Initial accomodation for iTPC
35 //
36 // Revision 2.99 2016/06/13 20:31:10 genevb
37 // Resolve Coverity BUFFER_SIZE_WARNING with careful copy function
38 //
39 // Revision 2.98 2016/06/10 02:55:54 genevb
40 // Coverity: memory leaks, possible null pointer dereferences, over-write character buffers
41 //
42 // Revision 2.97 2016/03/16 20:39:21 genevb
43 // remove accidental extraneous line
44 //
45 // Revision 2.96 2016/03/16 20:34:43 genevb
46 // Histogram list by subsystem, single TPC sector reference choice, and a couple histogram minima set
47 //
48 // Revision 2.95 2015/01/21 17:30:33 genevb
49 // Provide histogram normalization
50 //
51 // Revision 2.94 2014/04/10 17:58:40 genevb
52 // Fix for dE/dx slope plot
53 //
54 // Revision 2.93 2014/02/20 20:16:19 genevb
55 // Adjust dE/dx slope hist range, handle ROOT change for 2D polar plots
56 //
57 // Revision 2.92 2014/01/30 19:44:06 genevb
58 // Additional TPC histogram for monitoring gas contamination
59 //
60 // Revision 2.91 2013/03/12 03:42:55 genevb
61 // typo correction
62 //
63 // Revision 2.90 2013/03/12 03:41:18 genevb
64 // handle fms/fpd naming for now
65 //
66 // Revision 2.89 2013/03/12 03:05:43 genevb
67 // Add FMS/FPD histograms for Run 13+
68 //
69 // Revision 2.88 2012/05/01 18:37:19 genevb
70 // ZCol for Vtx XY distribution
71 //
72 // Revision 2.87 2012/03/07 02:04:10 genevb
73 // Differentiate event counts for Reference
74 //
75 // Revision 2.86 2012/03/03 01:29:06 genevb
76 // Output found/total vertices
77 //
78 // Revision 2.85 2012/01/31 22:14:53 genevb
79 // QA Shift Mode, optimized for AutoQA Browser
80 //
81 // Revision 2.84 2011/05/24 20:50:43 genevb
82 // Allow limited graphics file printing
83 //
84 // Revision 2.83 2011/05/13 21:12:58 genevb
85 // Hide original prefixes in title when plotting reference hists
86 //
87 // Revision 2.82 2011/03/15 21:05:25 genevb
88 // TPC hit phi sector labels
89 //
90 // Revision 2.81 2011/02/23 20:56:56 genevb
91 // Default to general histograms for references in absence of trig typed
92 //
93 // Revision 2.80 2011/02/23 18:46:29 genevb
94 // Fixed a bug introduced in version 2.71 when adding hists
95 //
96 // Revision 2.79 2011/02/19 02:43:39 genevb
97 // Fix those missing consts
98 //
99 // Revision 2.78 2011/02/19 02:22:18 genevb
100 // Allow for specification of histogram usage by the required detector sets
101 //
102 // Revision 2.77 2011/02/07 20:25:26 genevb
103 // Allow for limiting detectors
104 //
105 // Revision 2.76 2011/01/28 18:47:55 genevb
106 // Better handling of dirName
107 //
108 // Revision 2.75 2011/01/24 18:36:28 genevb
109 // Save hist list to results file even if not comparing
110 //
111 // Revision 2.74 2011/01/19 02:05:22 genevb
112 // Allow plain ROOT files with hists, and individual plot generation from 1 file
113 //
114 // Revision 2.73 2010/12/13 16:23:26 genevb
115 // Remove previous change, add Kolmogorov maximum distance as default
116 //
117 // Revision 2.72 2010/12/11 23:28:00 genevb
118 // Allow negative modes for 1-result
119 //
120 // Revision 2.71 2010/04/19 20:43:49 genevb
121 // Fixed bug with AddHists introduced in last fix
122 //
123 // Revision 2.70 2010/04/19 19:11:13 genevb
124 // Fixed bug with AddHists when some files are missing hists
125 //
126 // Revision 2.69 2010/04/09 21:13:28 genevb
127 // Use hobj pointer to ensure proper handling with reference hists
128 //
129 // Revision 2.68 2010/03/17 02:53:06 genevb
130 // Add hists even if not in first file
131 //
132 // Revision 2.67 2010/03/12 17:28:15 genevb
133 // Same ROOT quirk fix as done in rev. 2.64, but for reference hists
134 //
135 // Revision 2.66 2010/03/12 07:29:05 genevb
136 // Additional capability for saving images of each pad
137 //
138 // Revision 2.65 2010/03/08 18:04:33 genevb
139 // Include analysis score/result on plots
140 //
141 // Revision 2.64 2010/01/14 19:29:53 genevb
142 // Fix ROOT quirk with 1 page print, fix string/char conversions, protect LOG calls
143 //
144 // Revision 2.63 2009/05/11 17:52:14 genevb
145 // Fix fit-not-attached-to-hist-for-empty-hists, needed for ROOT 5.22
146 //
147 // Revision 2.62 2009/05/06 17:19:09 genevb
148 // Avoid PDF bug with text at angles
149 //
150 // Revision 2.61 2009/05/05 23:16:12 genevb
151 // Draw TPC sector boundaries and labels
152 //
153 // Revision 2.60 2009/05/05 00:31:07 genevb
154 // Add Anode guide lines in TPC Sector plots
155 //
156 // Revision 2.59 2009/05/04 23:37:40 genevb
157 // Add RDO boundary lines in TPC Sector plots
158 //
159 // Revision 2.58 2009/04/18 02:55:06 genevb
160 // Larger arrays for more trigger type hists
161 //
162 // Revision 2.57 2009/04/05 14:37:59 genevb
163 // Catch missing files
164 //
165 // Revision 2.56 2009/03/27 21:18:36 genevb
166 // Add Jet Patch trigger histograms
167 //
168 // Revision 2.55 2009/03/25 02:42:30 genevb
169 // Box --> Col for 2D plots
170 //
171 // Revision 2.54 2009/03/23 23:50:16 genevb
172 // Add color and logarithmic scales to TPC hit location plots
173 //
174 // Revision 2.53 2009/03/19 17:43:42 genevb
175 // Catch XYTpc hists which were polar
176 //
177 // Revision 2.52 2009/03/19 01:08:08 genevb
178 // Show both xy and rphi TPC hit hists
179 //
180 // Revision 2.51 2009/03/13 21:45:40 genevb
181 // Remove unhelpful stats
182 //
183 // Revision 2.50 2009/03/13 19:27:24 genevb
184 // Now draw TPC xy hits in polar coords
185 //
186 // Revision 2.49 2009/01/17 01:48:54 genevb
187 // Fixed broken multi-page output introduced in previous commit
188 //
189 // Revision 2.48 2009/01/08 23:40:13 genevb
190 // Introduce analyses with reference histograms
191 //
192 // Revision 2.47 2008/07/10 21:25:08 genevb
193 // Add canvas-to-code output with .CC suffix
194 //
195 // Revision 2.46 2008/07/09 20:52:16 genevb
196 // allow saving histograms to plain ROOT file
197 //
198 // Revision 2.45 2008/05/30 05:48:03 genevb
199 // Add ability to write out histogram data to .C files
200 //
201 // Revision 2.44 2008/05/28 05:16:06 genevb
202 // Allow summing over (ignoring) histogram prefixes
203 //
204 // Revision 2.43 2008/05/23 17:54:54 genevb
205 // Allow subset histogram list when copying/extracting
206 //
207 // Revision 2.42 2008/05/23 17:09:22 genevb
208 // Allow the use of any file for PrintList specification
209 //
210 // Revision 2.41 2008/05/15 19:37:17 genevb
211 // Changes to FTPC radial step plots/fits
212 //
213 // Revision 2.40 2008/02/06 23:31:32 genevb
214 // missing reset of logZ
215 //
216 // Revision 2.39 2007/12/14 02:18:29 genevb
217 // Add text of counts on NullPrimVtx plots
218 //
219 // Revision 2.38 2007/12/13 23:17:45 genevb
220 // Force 0 minimum on FTPC radial hists
221 //
222 // Revision 2.37 2007/04/28 20:36:15 perev
223 // Redundant StChain.h removed
224 //
225 // Revision 2.36 2007/04/24 00:33:06 genevb
226 // SSD hists; Logz for Dedx
227 //
228 // Revision 2.35 2007/04/20 03:42:35 genevb
229 // cout -> LOG_INFO
230 //
231 // Revision 2.34 2007/04/20 01:11:11 genevb
232 // ZCol on Dedx; printf -> LOG_INFO
233 //
234 // Revision 2.33 2007/04/07 04:39:04 genevb
235 // Use ZCol for PointXYTpc
236 //
237 // Revision 2.32 2007/03/13 18:42:27 genevb
238 // Add Svt list, simplified hlist include files, handle StMultiH2F, store dirName
239 //
240 // Revision 2.31 2007/02/26 20:45:00 genevb
241 // SVT drift hist
242 //
243 // Revision 2.30 2007/01/03 19:03:41 genevb
244 // Patch for hist titles removed after migration to Root vers. 5
245 //
246 // Revision 2.29 2006/05/18 16:38:03 genevb
247 // Introduce StHistUtil::GetRunYear()
248 //
249 // Revision 2.28 2006/04/24 21:23:13 genevb
250 // Fix problem with overlayed hists, and Patch for missing titles
251 //
252 // Revision 2.27 2006/03/28 21:35:31 genevb
253 // Single page output capability for eps,jpg,png,gif,tiff,etc. [see TPad::Print()]
254 //
255 // Revision 2.26 2006/03/28 01:58:38 genevb
256 // Allow PDF (and other) output formats (was only PostScript)
257 //
258 // Revision 2.25 2005/04/19 15:14:17 genevb
259 // Slight reordering of some FTPC code on user ranges in radial hists
260 //
261 // Revision 2.24 2005/02/08 17:12:37 genevb
262 // Limiting range on some PMD histos
263 //
264 // Revision 2.23 2004/12/13 15:52:35 genevb
265 // Numerous updates: PMD, primtrk, FPD, QAShift lists
266 //
267 // Revision 2.22 2004/10/04 16:40:41 genevb
268 // FTPC radial histos
269 //
270 // Revision 2.21 2004/06/09 22:01:31 genevb
271 // Modify line parameter for FTPC hist
272 //
273 // Revision 2.20 2004/05/24 15:13:10 genevb
274 // Range limit on FTPC radial step hists
275 //
276 // Revision 2.19 2004/02/12 16:54:22 genevb
277 // Separate MinBias histos
278 //
279 // Revision 2.18 2004/02/12 05:02:59 genevb
280 // Year 4 AuAu changes. New SVT histos.
281 //
282 // Revision 2.17 2004/02/10 16:31:15 genevb
283 // A few extra histo lines, features
284 //
285 // Revision 2.16 2003/09/02 17:55:26 perev
286 // gcc 3.2 updates + WarnOff
287 //
288 // Revision 2.15 2003/02/20 20:08:36 genevb
289 // Add new prefixes, other small modifications
290 //
291 // Revision 2.14 2003/02/15 22:00:14 genevb
292 // Add tpcSectors
293 //
294 // Revision 2.13 2003/01/22 21:32:46 genevb
295 // Fix Solaris compilation errors
296 //
297 // Revision 2.12 2002/09/06 02:51:34 genevb
298 // Remove limit on maximum number of histograms that can be copied
299 //
300 // Revision 2.11 2002/04/23 01:59:54 genevb
301 // Addition of BBC/FPD histos
302 //
303 // Revision 2.10 2002/02/12 18:41:57 genevb
304 // Additional FTPC histograms
305 //
306 // Revision 2.9 2002/01/26 03:04:05 genevb
307 // Fixed some problems with fcl histos
308 //
309 // Revision 2.8 2002/01/21 22:09:23 genevb
310 // Include some ftpc histograms from StFtpcClusterMaker
311 //
312 // Revision 2.7 2001/08/27 21:16:03 genevb
313 // Better filename hanlding
314 //
315 // Revision 2.6 2001/05/24 01:47:42 lansdell
316 // minor bug fixes and updates
317 //
318 // Revision 2.5 2001/05/23 00:14:27 lansdell
319 // added some logx code
320 //
321 // Revision 2.4 2001/04/25 14:17:16 genevb
322 // Reduced line width for newer Root
323 //
324 // Revision 2.3 2000/09/15 21:19:10 fisyak
325 // HPUX does not like delete [] newHist
326 //
327 // Revision 2.2 2000/08/28 19:21:40 genevb
328 // Plot MultiH1F hists like 1d
329 //
330 // Revision 2.1 2000/08/25 22:06:50 genevb
331 // Added histo descriptor in top right
332 //
333 // Revision 2.0 2000/08/25 15:47:38 genevb
334 // New revision: cleaned up, multiple PS files
335 //
336 //
338 // Histogram Utility methods for use with star makers and bfc output
340 
341 
342 #include <Stiostream.h>
343 #include <math.h>
344 #include <stdlib.h>
345 #include <string.h>
346 #include <Stsstream.h>
347 #include "TFile.h"
348 #include "TROOT.h"
349 #include "PhysicalConstants.h"
350 #include "TStyle.h"
351 #include "TCanvas.h"
352 #include "TObjString.h"
353 #include "TMath.h"
354 #include "TString.h"
355 #include "TRegexp.h"
356 #include "TPaveLabel.h"
357 #include "TPaveText.h"
358 #include "TLegend.h"
359 #include "TDatime.h"
360 #include "TLine.h"
361 #include "TLatex.h"
362 #include "StMessMgr.h"
363 
364 #include "St_DataSetIter.h"
365 #include "StMaker.h"
366 #include "TF1.h"
367 #include "TH3.h"
368 #include "TProfile.h"
369 
370 #include "StHistUtil.h"
371 
372 typedef TH1* TH1ptr;
373 typedef const char* charptr;
374 
375 const char* possibleQAPrefixes[10] = {"","LM","MM","HM","HP","XX","MB","CL","HT","JP"};
376 const char* possibleQASuffixes[10] = {
377  "General",
378  "Low Mult",
379  "Mid Mult",
380  "High Mult",
381  "High Pt",
382  "Other Physics",
383  "MinBias",
384  "Central",
385  "High Tower",
386  "Jet Patch"
387 };
388 
389 enum QAprintModes {QAprintSet,
390  QAprintSetRef,
391  QAprintIndiv,
392  QAprintIndivRef};
393 UInt_t QAU1 = 1;
394 
395 int sizeOfCharPtr = sizeof(Char_t*);
396 int sizeOfTH1Ptr = sizeof(TH1*);
397 
398 //_____________________________________________________________________________
399 // copied from StdEdxY2Maker vers. 1.83,
400 // avoids loading that library just for these
401 Double_t StdEdxY2Maker_gaus2(Double_t *x, Double_t *p) {
402  Double_t NormL = p[0];
403  Double_t mu = p[1];
404  Double_t muP = mu + p[4];
405  Double_t sigma = p[2];
406  Double_t sigmaP = TMath::Sqrt(sigma*sigma + 0.101741*0.101741);
407  Double_t phi = p[3];
408  Double_t frac = TMath::Sin(phi);
409  frac *= frac;
410  return TMath::Exp(NormL)*((1 - frac)*TMath::Gaus(x[0],mu ,sigma ,kTRUE) +
411  frac *TMath::Gaus(x[0],muP,sigmaP,kTRUE));
412 }
413 TF1 *StdEdxY2Maker_Gaus2() {
414  TF1 *f = new TF1("Gaus2",StdEdxY2Maker_gaus2,-3,3,5);
415  f->SetParName(0,"NormL"); f->SetParLimits(0,-10.,10.);
416  f->SetParName(1,"mu"); f->SetParLimits(1,-0.5,0.5);
417  f->SetParName(2,"sigma"); f->SetParLimits(2, 0.2,0.5);
418  f->SetParName(3,"phiP"); f->SetParLimits(3, 0.0,TMath::Pi()/4);
419  f->SetParName(4,"muP");
420  f->SetParameters(10,0,0.3,0.1,1.315);
421  // f->FixParameter(4,1.425);
422  return f;
423 }
424 //_____________________________________________________________________________
425 
426 ClassImp(StHistUtil)
427 
428 //_____________________________________________________________________________
429 
430 // Constructor
431 
433 
434  numOfPosPrefixes = 10;
435  possiblePrefixes = possibleQAPrefixes;
436  possibleSuffixes = possibleQASuffixes;
437 
438  m_ListOfLogY = 0;
439  m_ListOfLogX = 0;
440  m_ListOfPrint = 0;
441  m_HistCanvas = 0;
442  m_HistCanvasR = 0;
443  debug = kFALSE;
444  m_CurPrefix = -1;
445  m_CurPage = 0;
446  m_CurFileName = "";
447  m_OutFileName = "";
448  m_OutType = "ps"; // postscript output by default
449  m_PrintMode = 0;
450  m_OutMultiPage = kTRUE;
451  m_OutIndividuals = "";
452  m_QAShiftMode = kFALSE;
453  m_RunYear = 0;
454 
455  Ltitle = 0;
456  Ldesc = 0;
457 
458  maxHistCopy = 4096;
459  newHist = new TH1ptr[maxHistCopy];
460  memset(newHist,0,maxHistCopy*sizeOfTH1Ptr);
461  m_dirName[0] = 0;
462 
463  ignorePrefixes = kFALSE;
464 
465  m_analMode = kFALSE;
466  m_refResultsFile[0] = 0;
467  m_refOutFile[0] = 0;
468  m_refCuts = 0;
469  m_refInFile = 0;
470  m_PntrToMaker = 0;
471  m_PntrToPlainFile = 0;
472 
473  m_PadColumns = 0;
474  m_PadRows = 0;
475  m_PaperWidth = 0;
476  m_PaperHeight = 0;
477 
478  m_Detectors = "";
479 
480  m_ItemsToClear = new TList();
481  m_ItemsToClear->SetOwner();
482 }
483 //_____________________________________________________________________________
484 
485 // Destructor
486 
487 StHistUtil::~StHistUtil(){
488  SafeDelete(m_HistCanvas);
489  SafeDelete(m_HistCanvasR);
490  if (m_ListOfLogY) {
491  m_ListOfLogY->Delete();
492  SafeDelete(m_ListOfLogY);
493  }
494  if (m_ListOfLogX) {
495  m_ListOfLogX->Delete();
496  SafeDelete(m_ListOfLogX);
497  }
498  if (m_ListOfPrint) {
499  m_ListOfPrint->Delete();
500  SafeDelete(m_ListOfPrint);
501  }
502  if (newHist) {
503  for (int ijk=0; ijk<maxHistCopy; ijk++) delete newHist[ijk];
504  delete [] newHist;
505  }
506  delete m_ItemsToClear;
507 }
508 //_____________________________________________________________________________
509 void StHistUtil::Clear() {
510  m_ItemsToClear->Clear();
511 }
512 //_____________________________________________________________________________
513 void StHistUtil::SetOutFile(const Char_t *fileName, const Char_t* type) {
514  m_OutFileName = fileName;
515  if (m_OutFileName.EndsWith("+")) {
516  m_OutIndividuals = ".eps"; // only option working in ROOT 5.22
517  m_PrintMode |= QAU1<<QAprintIndiv;
518  if (m_OutFileName.EndsWith("++")) m_PrintMode |= QAU1<<QAprintIndivRef;
519  while (m_OutFileName.EndsWith("+")) m_OutFileName.Chop();
520  }
521  if (type) {
522  m_OutType = type;
523  } else {
524  if (m_OutFileName.EndsWith(".ps")) m_OutType="ps";
525  else if (m_OutFileName.EndsWith(".eps")) m_OutType="eps";
526  else if (m_OutFileName.EndsWith(".epsf")) m_OutType="Preview";
527  else if (m_OutFileName.EndsWith(".pdf")) m_OutType="pdf";
528  else if (m_OutFileName.EndsWith(".jpg")) m_OutType="jpg";
529  else if (m_OutFileName.EndsWith(".jpeg")) m_OutType="jpg";
530  else if (m_OutFileName.EndsWith(".gif")) m_OutType="gif";
531  else if (m_OutFileName.EndsWith(".tif")) m_OutType="tiff";
532  else if (m_OutFileName.EndsWith(".tiff")) m_OutType="tiff";
533  else if (m_OutFileName.EndsWith(".svg")) m_OutType="svg";
534  else if (m_OutFileName.EndsWith(".xpm")) m_OutType="xpm";
535  else if (m_OutFileName.EndsWith(".png")) m_OutType="png";
536  else if (m_OutFileName.EndsWith(".CC")) m_OutType="CC"; // Save canvases as code
537  else if (m_OutFileName.EndsWith(".C")) m_OutType="C"; // Save histograms as code
538  else if (m_OutFileName.EndsWith(".root")) m_OutType="root";
539  else if (m_OutFileName.EndsWith("none")) m_OutType="none"; // No output set
540  else if (m_OutFileName.EndsWith(".qas")) m_OutType="qas"; // QA Shift mode
541  else {
542  LOG_INFO << "SetHistUtil::SetOutFile(): unknown type, assuming ps" << endm;
543  m_OutType = "ps";
544  m_OutFileName.Append(".ps");
545  }
546  }
547  if (!m_OutType.CompareTo("qas")) { // QA Shift mode options
548  m_QAShiftMode = kTRUE;
549  m_OutIndividuals = ".svg";
550  m_PrintMode |= QAU1<<QAprintIndiv;
551  m_PrintMode |= QAU1<<QAprintIndivRef;
552  m_OutType = "none";
553  }
554  if (m_OutType.CompareTo("none")) {
555  m_PrintMode |= QAU1<<QAprintSet;
556  m_PrintMode |= QAU1<<QAprintSetRef;
557  }
558 
559  // Multipage output for ps,pdf
560  m_OutMultiPage = !(m_OutType.CompareTo("ps")
561  && m_OutType.CompareTo("pdf") );
562  if (m_OutMultiPage) {
563  LOG_INFO << "StHistUtil::SetOutFile(): Multipage output" << endm;
564  } else {
565  LOG_INFO << "StHistUtil::SetOutFile(): Single page output" << endm;
566  }
567 }
568 //_____________________________________________________________________________
569 void StHistUtil::CloseOutFile() {
570  m_HistCanvas->Modified();
571  m_HistCanvas->Update();
572  if (!m_CurFileName.IsNull()) {
573  if (m_OutMultiPage) m_CurFileName.Append(")");
574  if (m_OutType.CompareTo("CC")) {
575  // single page seems to have a bug with "()" notation as of Root 5.22.00
576  if (m_CurPage==1) m_CurFileName.Chop().Chop();
577  if (m_PrintMode & QAU1<<QAprintSet)
578  m_HistCanvas->Print(m_CurFileName.Data(),m_OutType.Data());
579  } else
580  m_HistCanvas->SaveSource(m_CurFileName.Data());
581  if (m_refInFile) {
582  m_HistCanvasR->Modified();
583  m_HistCanvasR->Update();
584  m_CurFileNameR.Append(")");
585  // single page seems to have a bug with "()" notation as of Root 5.22.00
586  if (m_CurPage==1) m_CurFileNameR.Chop().Chop();
587  if (m_PrintMode & QAU1<<QAprintSetRef)
588  m_HistCanvasR->Print(m_CurFileNameR.Data(),m_OutType.Data());
589  // anal mode doesn't support single page output
590  }
591  } else {
592  LOG_INFO << "StHistUtil::CloseOutFile(): No output file" << endm;
593  }
594 }
595 //_____________________________________________________________________________
596 TString StHistUtil::StripPrefixes(const Char_t* histName, Int_t& prenum, Int_t mode) {
597  // mode < 0 : strip maker name
598  // mode > 0 : strip trigger prefix
599  // mode = 0 : strip both
600  // Figure out and strip appropriate prefix index
601  TString hName(histName);
602  Char_t makerBuffer[4];
603  memset(makerBuffer,0,4);
604  if ((hName.BeginsWith("Tab")) || (hName.BeginsWith("StE"))) {
605  memcpy(makerBuffer,histName,3);
606  hName.Remove(0,3);
607  }
608  prenum = 0; // Possible prefix=0 means no prefix
609  if (mode >= 0) {
610  for (Int_t i=1; i<numOfPosPrefixes; i++) {
611  if (hName.BeginsWith(possiblePrefixes[i])) {
612  prenum = i;
613  hName.Remove(0,strlen(possiblePrefixes[i]));
614  break;
615  }
616  }
617  if (mode>0) hName.Prepend(makerBuffer);
618  }
619  return hName;
620 }
621 //_____________________________________________________________________________
622 Bool_t StHistUtil::CheckOutFile(const Char_t *histName) {
623 // Method to determine appropriate PostScript file for output
624 
625  // Figure out appropriate prefix index
626  Int_t newPrefix = -1;
627  StripPrefixes(histName,newPrefix);
628 
629  if (newPrefix == m_CurPrefix) return kFALSE;
630 
631  CloseOutFile();
632  m_CurPrefix = newPrefix;
633  if (m_OutType.CompareTo("C") && m_OutType.CompareTo("root")) {
634  m_CurFileName = m_OutFileName;
635  Ssiz_t insertPos = m_CurFileName.Last('.');
636  if (insertPos<0) insertPos = m_CurFileName.Length();
637  if (m_OutMultiPage) m_CurFileName.Append("(");
638  else m_CurFileName.Insert(insertPos,"_");
639  m_CurFileName.Insert(insertPos,possiblePrefixes[m_CurPrefix]);
640  }
641  (m_CurFileNameR = "Ref_") += m_CurFileName;
642 
643  if (Ldesc) {
644  Ldesc->Clear();
645  Ldesc->AddText(possibleSuffixes[m_CurPrefix]);
646  Ldesc->AddText("Hists");
647  }
648  return kTRUE;
649 }
650 //_____________________________________________________________________________
651 Int_t StHistUtil::DrawHists(const Char_t *dirName) {
652 // Method DrawHists -->
653 // Plot the selected histograms and generate the postscript file as well
654 
655  LOG_INFO << " **** Now in StHistUtil::DrawHists **** " << endm;
656 
657  Int_t canvasWidth,canvasHeight;
658 
659  if (m_QAShiftMode) {
660  LOG_INFO << "In QA Shift Mode - overriding other inputs" << endm;
661  m_PadColumns=1;
662  m_PadRows=1;
663  canvasWidth = 250;
664  canvasHeight = 250;
665  } else {
666  // SetPaperSize wants width & height in cm: A4 is 20,26 & US is 20,24
667  gStyle->SetPaperSize(m_PaperWidth,m_PaperHeight);
668  // TCanvas wants width & height in pixels (712 x 950 corresponds to A4 paper)
669  // (600 x 780 US )
670  // TCanvas *m_HistCanvas = new TCanvas("CanvasName","Canvas Title",30*m_PaperWidth,30*m_PaperHeight);
671  canvasWidth = 600;
672  canvasHeight = 780;
673  }
674 
675  //set Style of Plots
676  const Int_t numPads = m_PadColumns*m_PadRows;
677  gStyle->SetStatStyle(0);
678  gStyle->SetOptDate(0);
679  gStyle->SetPalette(1);
680 
681 
682  //setup canvas
683  SafeDelete(m_HistCanvas);
684  SafeDelete(m_HistCanvasR);
685 
686  if (m_refInFile) {
687  m_HistCanvasR = new TCanvas("CanvasNameR"," STAR Reference Histogram Canvas",20,20,canvasWidth,canvasHeight);
688  }
689  m_HistCanvas = new TCanvas("CanvasName"," STAR Maker Histogram Canvas",0,0,canvasWidth,canvasHeight);
690 
691 
692  TPad *graphPad = m_HistCanvas;
693  TPad *graphPadR = m_HistCanvasR;
694  TPaveLabel* LtitleR = 0;
695  TPaveLabel *Ldatetime = 0;
696  TPaveLabel *Lpage = 0;
697  m_CurPage=1;
698 
699  if (!m_QAShiftMode) {
700  // Do not draw page numbers, titles, etc. in QA Shift mode
701 
702  // write title at top of canvas - first page
703  Ltitle = new TPaveLabel(0.08,0.96,0.88,1.0,m_GlobalTitle.Data(),"br");
704  Ltitle->SetFillColor(18);
705  Ltitle->SetTextFont(32);
706  Ltitle->SetTextSize(0.5);
707  Ltitle->Draw();
708  if (m_refInFile) {
709  m_HistCanvasR->cd();
710  LtitleR = new TPaveLabel(0.08,0.96,0.88,1.0,m_refInFile->GetName(),"br");
711  LtitleR->SetFillColor(18);
712  LtitleR->SetTextFont(32);
713  LtitleR->SetTextSize(0.5);
714  LtitleR->Draw();
715  m_HistCanvas->cd();
716  }
717 
718  // write descriptor at top of canvas - first page
719  Ldesc = new TPaveText(0.90,0.96,0.99,1.0,"br");
720  Ldesc->SetFillColor(18);
721  Ldesc->SetTextFont(32);
722  Ldesc->Draw();
723  if (m_refInFile) {
724  m_HistCanvasR->cd();
725  Ldesc->Draw();
726  m_HistCanvas->cd();
727  }
728 
729  // now put in date & time at bottom right of canvas - first page
730  TDatime HistTime;
731  const Char_t *myTime = HistTime.AsString();
732  Ldatetime = new TPaveLabel(0.7,0.01,0.95,0.03,myTime,"br");
733  Ldatetime->SetTextSize(0.6);
734  Ldatetime->Draw();
735  if (m_refInFile) {
736  m_HistCanvasR->cd();
737  Ldatetime->Draw();
738  m_HistCanvas->cd();
739  }
740 
741 
742  // now put in page # at bottom left of canvas - first page
743  Lpage = new TPaveLabel(0.1,0.01,0.16,0.03,Form("%d",m_CurPage),"br");
744  Lpage->SetTextSize(0.6);
745  Lpage->Draw();
746  if (m_refInFile) {
747  m_HistCanvasR->cd();
748  Lpage->Draw();
749  m_HistCanvas->cd();
750  }
751 
752  // Make 1 big pad on the canvas - make it a little bit inside the canvas
753  // - must cd to get to this pad!
754  // order is x1 y1 x2 y2
755  graphPad = new TPad("PadName","Pad Title",0.0,0.05,1.00,0.95);
756  graphPad->Draw();
757  graphPad->cd();
758  graphPad->Divide(m_PadColumns,m_PadRows);
759  if (m_refInFile) {
760  m_HistCanvasR->cd();
761  graphPadR = new TPad("PadNameR","Pad TitleR",0.0,0.05,1.00,0.95);
762  graphPadR->Draw();
763  graphPadR->cd();
764  graphPadR->Divide(m_PadColumns,m_PadRows);
765  graphPad->cd();
766  }
767 
768  } // !m_QAShiftMode
769 
770  Int_t padCount = 0;
771  Bool_t padAdvance = kTRUE;
772 
773 
774  // Now find the histograms
775  // get the TList pointer to the histograms:
776  PathCopy(m_dirName,dirName);
777  TList* dirList = (m_PntrToMaker ? FindHists(m_dirName) : FindHists(m_PntrToPlainFile));
778  if (!dirList) { LOG_INFO << " DrawHists - histograms not available! " << endm; }
779 
780  TIter nextHist(dirList);
781  Int_t histCounter = 0;
782  Int_t histReadCounter = 0;
783  Bool_t started = kFALSE;
784 
785 
786  TObject *obj = 0;
787  TLine ruler;
788  TLatex latex;
789  TF1* fcsFunc1 = 0;
790  TF1* fcsFunc2 = 0;
791  TF1* fcsFunc3 = 0;
792 
793  ofstream* C_ostr = 0;
794  TFile* root_ofile = 0;
795  if (!m_OutType.CompareTo("C")) {
796  C_ostr = new ofstream(m_OutFileName);
797  (*C_ostr) << " gSystem->Load(\"St_base\");" << endl;
798  (*C_ostr) << " gSystem->Load(\"StUtilities\");" << endl;
799  // not supporting ref output
800  } else if (!m_OutType.CompareTo("root"))
801  root_ofile = new TFile(m_OutFileName,"RECREATE");
802 
803  // Reference analyses:
804  TList* dirListR = 0;
805  ofstream* R_ostr = 0;
806  if (m_analMode) {
807  dirListR = FindHists(m_refInFile);
808  // By default, save histograms to a ROOT file as future reference
809  if (!root_ofile) root_ofile = new TFile(m_refOutFile,"RECREATE");
810  }
811 
812  // function to fit FTPC radial step
813  static TF1* fitFRS = 0;
814  if (!fitFRS) fitFRS = new TF1("fitFRS","[0]*(x<[1]-[2])+([0]+([3]-[0])*(x-([1]-[2]))/(2*[2]))*(x>[1]-[2])*(x<[1]+[2])+[3]*(x>[1]+[2])",6.5, 9.0);
815 
816  while ((obj = nextHist())) {
817 
818  if (!obj->InheritsFrom("TH1")) continue;
819 
820  TH1* hobj = (TH1*) obj;
821  const char* oname = hobj->GetName();
822  const char* otitle = hobj->GetTitle();
823  TString oName = oname;
824  oName.ReplaceAll(' ','_');
825  histReadCounter++;
826  LOG_INFO << Form(" %d. Reading ... %s::%s; Title=\"%s\"\n",
827  histReadCounter,hobj->ClassName(),oname, otitle) << endm;
828  if (!started && (m_FirstHistName.CompareTo("*")==0 ||
829  m_FirstHistName.CompareTo(oName)==0))
830  started = kTRUE;
831 
832 // Here is the actual loop over histograms !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
833  if (started) {
834  if (oName.CompareTo(m_LastHistName)==0) started = kFALSE;
835  histCounter++;
836 //...........................................................................
837 
838 
839  // If there's no print list, then print all histograms in directory
840  // If there is a print list, then only print if hist name is on list
841  if (!m_ListOfPrint || (m_ListOfPrint->FindObject(oname))) {
842 
843  // Some ROOT output options
844  if (root_ofile) {
845  root_ofile->cd();
846  hobj->Write();
847  }
848  if (C_ostr) {
849  hobj->SavePrimitive(*C_ostr);
850  continue;
851  }
852 
853  // this histogram will actually be printed/drawn!!
854  LOG_INFO << Form(" - %d. Drawing ... %s::%s; Title=\"%s\"\n",
855  histCounter,hobj->ClassName(),oname, otitle) << endm;
856 
857  // Switch to a new page...............................
858  if (CheckOutFile(oname)) {
859  padCount = numPads;
860  m_CurPage = 0;
861  }
862  if (m_QAShiftMode) {
863  graphPad->Clear();
864  if (m_refInFile) graphPadR->Clear();
865  padCount = 0;
866  m_CurPage++;
867  } else if (padCount == numPads) {
868  // must redraw the histcanvas for each new page!
869  m_HistCanvas->Modified();
870  m_HistCanvas->Update();
871  if (m_PrintMode & QAU1<<QAprintSet &&
872  m_CurPage>0 && !m_CurFileName.IsNull()) {
873  if (m_OutType.CompareTo("CC")) {
874  m_HistCanvas->Print(m_CurFileName.Data(),m_OutType.Data());
875  } else
876  m_HistCanvas->SaveSource(m_CurFileName.Data());
877  m_CurFileName.ReplaceAll("(",0); // doesn't hurt to do > once
878  } else {
879  m_HistCanvas->Draw();
880  }
881 
882  while (padCount > 0) graphPad->GetPad(padCount--)->Clear();
883 
884  if (m_refInFile) {
885  m_HistCanvasR->Modified();
886  m_HistCanvasR->Update();
887  if (m_PrintMode & QAU1<<QAprintSetRef &&
888  m_CurPage>0 && !m_CurFileName.IsNull()) {
889  m_HistCanvasR->Print(m_CurFileNameR.Data(),m_OutType.Data());
890  m_CurFileNameR.ReplaceAll("(",0); // doesn't hurt to do > once
891  } else m_HistCanvasR->Draw();
892  padCount = numPads;
893  while (padCount > 0) graphPadR->GetPad(padCount--)->Clear();
894  }
895 
896  // update the page number
897  m_CurPage++;
898  if (Lpage) Lpage->SetLabel(Form("%d",m_CurPage));
899 
900  if (!m_OutMultiPage && !m_CurFileName.IsNull()) {
901  Ssiz_t last_us = m_CurFileName.Last('_') + 1;
902  Ssiz_t lastdot = m_CurFileName.Last('.') - last_us;
903  m_CurFileName.Replace(last_us,lastdot,Form("%d",m_CurPage));
904  }
905  }
906 
907  // check dimension of histogram
908  Int_t chkdim = hobj->GetDimension();
909 
910  // go to next pad
911  int curPad = (++padCount);
912 
913  // Will need to display both histograms if doing reference analysis...
914  TH1* hobjO = hobj;
915  TH1* hobjR = 0;
916  if (dirListR) {
917  Int_t tempint = -1;
918  // try: full name
919  TString onamebase = oname;
920 #define SingleTpcSectorReference false
921  if (SingleTpcSectorReference && onamebase.Contains("iTpcSector")) {
922  // last parameter is the single sector to use for reference:
923  // e.g. TpcSector14 => TpcSector20 if the number is "20"
924  onamebase.Replace(onamebase.Index("iTpcSector")+10,2,"16");
925  } else if (SingleTpcSectorReference && onamebase.Contains("TpcSector")) {
926  // last parameter is the single sector to use for reference:
927  // e.g. TpcSector20 => TpcSector14 if the number is "14"
928  onamebase.Replace(onamebase.Index("TpcSector")+9,2,"9");
929  }
930  hobjR = (TH1*) (dirListR->FindObject(onamebase.Data()));
931  if (!hobjR) {
932  // try: strip just maker from name
933  onamebase = StripPrefixes(oname,tempint,-1);
934  hobjR = (TH1*) (dirListR->FindObject(onamebase.Data()));
935  if (!hobjR) {
936  // try: strip just trigger type from name
937  onamebase = StripPrefixes(oname,tempint,1);
938  hobjR = (TH1*) (dirListR->FindObject(onamebase.Data()));
939  if (!hobjR) {
940  // try: strip maker and trigger type from name
941  onamebase = StripPrefixes(oname,tempint,0);
942  hobjR = (TH1*) (dirListR->FindObject(onamebase.Data()));
943  }
944  }
945  }
946  if (hobjR) {
947  TString htitle = StripPrefixes(hobjR->GetTitle(),tempint,0);
948  if (!htitle.BeginsWith("Ref:")) htitle.Prepend("Ref:");
949  hobjR->SetTitle(htitle.Data());
950  }
951  }
952  TVirtualPad* objPad = 0;
953  for (int analRepeat = 0;analRepeat < (hobjR ? 2 : 1); analRepeat++) {
954 
955  padAdvance = kTRUE;
956  if (analRepeat) {
957  objPad=gPad;
958  graphPadR->cd(m_QAShiftMode ? 0 : curPad);
959  hobj=hobjR;
960  } else graphPad->cd(m_QAShiftMode ? 0 : curPad);
961 
962  // set x & y grid off by default
963  gPad->SetGridy(0);
964  if (oName.Contains("H_matchCand")) {
965  gPad->SetGridx(1);
966  gStyle->SetGridStyle(6);
967  gStyle->SetGridColor(kOrange);
968  } else {
969  gPad->SetGridx(0);
970  gStyle->SetGridStyle(3);
971  gStyle->SetGridColor(kGray);
972  }
973 
974  // set stats to draw
975  if (oName.Contains("TpcSector") ||
976  oName.Contains("PointRPTpc") ||
977  oName.Contains("PointXYTpc") ||
978  oName.Contains("TrigBits")) {
979  gStyle->SetOptStat(11);
980  } else if (oName.Contains("NullPrim")) {
981  gStyle->SetOptStat(1111);
982  } else {
983  gStyle->SetOptStat(111111);
984  }
985 
986  // set bottom margin of pad
987  if (oName.Contains("TrigBits")) {
988  gPad->SetBottomMargin(0.35);
989  hobj->GetXaxis()->SetLabelSize(0.03);
990  hobj->GetXaxis()->SetLabelFont(42);
991  } else {
992  gPad->SetBottomMargin(0.10);
993  }
994 
995  if (oName.Contains("GtrkPadfT")) hobj->SetMinimum(0.8);
996 
997  // set logX,Y,Z scale on/off
998 
999 // Normalized histograms: rely on negative bin content 0 as a flag (seems OK for now, but not guaranteed)
1000  Float_t BinCont0 = hobj->GetBinContent(0);
1001  if (BinCont0 < 0) {
1002  LOG_INFO << " -- Will normalize by 1/" << -BinCont0 << ": " << oname << endm;
1003  hobj->Scale(-1.0/BinCont0);
1004  }
1005 
1006 
1007 // Set logY scale on if: there is a loglist, if the hist name is on the list, if it has entries
1008 // and if the max entries in all bins is > 0
1009  if (m_ListOfLogY && m_ListOfLogY->FindObject(oname) &&
1010  hobj->GetEntries() && hobj->GetMaximum() ) {
1011  gPad->SetLogy(1);
1012  if (!analRepeat) {LOG_INFO << " -- Will draw in logY scale: " << oname <<endm;}
1013  } else gPad->SetLogy(0);
1014 
1015 
1016 // Set logX scale on if: there is a loglist, if the hist name is on the list, if it has entries
1017 // and if the max entries in all bins is > 0
1018  if (m_ListOfLogX && m_ListOfLogX->FindObject(oname) &&
1019  hobj->GetEntries() && hobj->GetMaximum() ) {
1020  gPad->SetLogx(1);
1021  if (!analRepeat) {LOG_INFO << " -- Will draw in logX scale: " << oname <<endm;}
1022  } else gPad->SetLogx(0);
1023 
1024 // Set logZ scale
1025  if (oName.EndsWith("PVsDedx") ||
1026  oName.Contains("fms_qt_") ||
1027  oName.Contains("fpd_channel_") ||
1028  oName.Contains("TofPID") ||
1029  oName.Contains("multiplicity_etofHits") ||
1030  oName.Contains("eTofHits") ||
1031  oName.Contains("etofMult") ||
1032  oName.Contains("G_matchCand_") ||
1033  oName.Contains("RP_cluster_xy") ||
1034  oName.Contains("TpcSector") ||
1035  oName.Contains("PointRPTpc") ||
1036  oName.Contains("PointXYTpc")) {
1037  gPad->SetLogz(1);
1038  if (!analRepeat) {LOG_INFO << " -- Will draw in logZ scale: " << oname <<endm;}
1039  } else gPad->SetLogz(0);
1040 
1041 // Limit x range for some histograms
1042  if (oName.EndsWith("QaPointTpc") ||
1043  oName.EndsWith("QaPointSvt") ||
1044  oName.Contains("QaPointSvtLaser") ||
1045  oName.EndsWith("QaPointSsd") ||
1046  oName.EndsWith("QaPointFtpc") ||
1047  oName.EndsWith("QaRichTot") ||
1048  oName.EndsWith("QaV0Vtx") ||
1049  oName.EndsWith("QaXiVtxTot") ||
1050  oName.Contains("QaPmdTotal") ||
1051  oName.Contains("QaCpvTotal") ||
1052  oName.EndsWith("SImpactTime") ||
1053  oName.EndsWith("trkGoodTTS")) {
1054  Float_t mean = hobj->GetMean(1);
1055  Float_t window = hobj->GetRMS(1);
1056  Float_t bwid = hobj->GetXaxis()->GetBinWidth(1);
1057  if (window < bwid) window = bwid;
1058  hobj->SetAxisRange(mean-5*window,mean+5*window,"X");
1059  }
1060 
1061 // Limit both x & y ranges together for some histograms
1062  if (oName.EndsWith("trkGoodF")) {
1063  Float_t mean1 = hobj->GetMean(1);
1064  Float_t mean2 = hobj->GetMean(2);
1065  Float_t window1 = hobj->GetRMS(1);
1066  Float_t window2 = hobj->GetRMS(2);
1067  Float_t bwid = hobj->GetXaxis()->GetBinWidth(1);
1068  if (window1 < bwid) window1 = bwid;
1069  if (window2 < bwid) window2 = bwid;
1070  Float_t lo = TMath::Min(mean1-5*window1,mean2-5*window2);
1071  Float_t hi = TMath::Max(mean1+5*window1,mean2+5*window2);
1072  hobj->SetAxisRange(lo,hi,"X");
1073  hobj->SetAxisRange(lo,hi,"Y");
1074  }
1075 
1076 // Limit both x & y ranges for QaVtxFtpcE/WTpcXY histograms for non-year 7 runs
1077  if (oName.Contains("VtxFtpc")&&oName.Contains("TpcXY")&&m_RunYear!=7){
1078  hobj->GetXaxis()->SetRangeUser(-2.0,2.0);
1079  hobj->GetYaxis()->SetRangeUser(-2.0,2.0);
1080  }
1081 
1082 // Limit x (time) range for some ETOF histograms
1083  if (oName.Contains("Diff_etofHits")) {
1084  hobj->GetXaxis()->SetRangeUser(0.0,100.0);
1085  }
1086 
1087  // actually draw,print
1088  if ((chkdim == 3) && (hobj->InheritsFrom("StMultiH2F"))) {
1089  hobj->Draw("Col");
1090  } else if ((chkdim == 3) && (oName.Contains("Z3A"))) {
1091  latex.SetTextAngle(0);
1092  latex.SetTextAlign(12);
1093  TH3F* Z3A = (TH3F*) hobj;
1094  Bool_t noneYet = kTRUE;
1095  // copied from StdEdxY2Maker::FinishRun() vers. 1.83,
1096  // with minor modifications for plotting the results
1097  Double_t slope = 1.7502e-6;// slope from Blair 1/( O2 in ppm., cm )
1098  const Char_t *IO[2] = {"Inner", "Outer"};
1099  Double_t xmin[2] = { 40, 40};
1100  Double_t xmax[2] = {200,180};
1101  TF1 *gg = StdEdxY2Maker_Gaus2();
1102  float histmiddle = 0;
1103  for (Int_t io = 1; io <= 2; io++) {
1104  Z3A->GetXaxis()->SetRange(io,io);
1105  TH2 *I = (TH2 *) Z3A->Project3D(Form("zy%i",io));
1106  if (I) {
1107  I->FitSlicesY(gg);
1108  TH1D *proj = (TH1D*) gDirectory->Get(Form("%s_1",I->GetName()));
1109  if (proj) {
1110  proj->Fit("pol1","erq",(noneYet ? "" : "same"),xmin[io-1],xmax[io-1]);
1111  proj->SetLineColor(8-io);
1112  proj->SetMarkerColor(8-io);
1113  proj->SetMarkerStyle(24-io);
1114  proj->SetStats(0);
1115  proj->SetTitle(otitle);
1116  TF1 *f = (TF1 *) proj->GetListOfFunctions()->FindObject("pol1");
1117  if (f) {
1118  gMessMgr->Info() << "StHistUtil: Estimated content of O2 (ppm) "
1119  << "from slope in drift distance for "
1120  << Form("%s = %10.2f +/- %10.2f", IO[io-1],
1121  -f->GetParameter(1)/slope, f->GetParError(1)/slope)
1122  << endm;
1123  if (noneYet) histmiddle = f->Eval(100.0);
1124  f->SetLineColor(6-2*io);
1125  latex.SetTextColor(6-2*io);
1126  latex.DrawLatex(20,histmiddle+0.9-0.6*io,Form("%s : %10.2f +/- %10.2f\n",
1127  IO[io-1],-f->GetParameter(1)/slope,f->GetParError(1)/slope));
1128  }
1129  if (noneYet) {
1130  proj->SetMinimum(histmiddle-0.4);
1131  proj->SetMaximum(histmiddle+0.4);
1132  noneYet = kFALSE;
1133  }
1134  }
1135  }
1136  }
1137  latex.SetTextColor(1);
1138  } else if ((chkdim==2) && (oName.Contains("PointRPTpc") ||
1139  (oName.Contains("PointXYTpc") && // Unfortunately was polar for a short time
1140  TMath::Abs((hobj->GetYaxis()->GetXmax()/TMath::Pi())-2)<1e-5))) {
1141  TH2F* htmp = new TH2F(Form("%s.",hobj->GetName()),hobj->GetTitle(),1,-200,200,1,-200,200);
1142  float hmin = (oName.Contains("PointRPTpcQ") ? 1e-4 : 1.0);
1143  htmp->Fill(0.,0.,hmin); htmp->SetMinimum(hmin);
1144  htmp->SetStats(kFALSE);
1145  htmp->Draw();
1146  hobj->SetMinimum(0.9*hmin);
1147  if (gROOT->GetVersionInt() < 52800) {
1148  hobj->Draw("Pol ZCol Same");
1149  } else {
1150  // lego plots always needed phi,r from x,y
1151  // (z)col plots, however, needed r,phi from x,y
1152  // Now, (z)col plots also need phi,r from x,y
1153  // https://sft.its.cern.ch/jira/browse/ROOT-2845
1154  FlipAxes(hobj)->Draw("Pol ZCol Same");
1155  }
1156  } else if ((chkdim == 2) &&
1157  (oName.EndsWith("SImpactTime"))) {
1158  hobj->SetMarkerStyle(7);
1159  ((TH2F*) hobj)->Rebin2D(100,1,0);
1160  hobj->SetXTitle("Time in run [sec]");
1161  hobj->SetYTitle("signed impact (sDCA) [cm]");
1162  ((TH2F*) hobj)->ProfileX()->Draw();
1163  } else if ((chkdim == 2) &&
1164  (oName.EndsWith("SvtLoc") ||
1165  oName.EndsWith("PVsDedx") ||
1166  oName.EndsWith("VtxPrXY") ||
1167  oName.EndsWith("SSD") ||
1168  oName.EndsWith("PointXYSvt") ||
1169  oName.Contains("TpcSector") ||
1170  oName.Contains("PointXYTpc") ||
1171  oName.Contains("SectorvsSensor") ||
1172  oName.Contains("LaddervsSensor"))) {
1173  if (oName.Contains("TpcSector")) {
1174  // method to suppress a single hot channel
1175  Double_t max1 = hobj->GetMaximum();
1176  Double_t max2 = hobj->GetMaximum(max1);
1177  if (max1/max2 > 5) hobj->SetBinContent(hobj->GetMaximumBin(),max2);
1178  }
1179  hobj->Draw("ZCol");
1180  } else if ((chkdim == 2) && (!hobj->InheritsFrom("StMultiH1F"))) {
1181  hobj->Draw("Col");
1182  if ((oName.EndsWith("trkGoodF"))||(oName.EndsWith("VtxSvtvsTpc"))) {
1183  ruler.SetLineColor(46);
1184  ruler.SetLineWidth(2);
1185  ruler.DrawLineNDC(0.1,0.1,0.9,0.9);
1186  }
1187  } else {
1188  if (oName.Contains("QaBbc") ||
1189  (oName.Contains("QaPmd") && !oName.Contains("Total")) ||
1190  (oName.Contains("QaFpd") && !oName.Contains("Sums"))) {
1191  hobj->SetBarOffset();
1192  }
1193  if (oName.Contains("TrigBits")) {
1194  hobj->SetLineWidth(1);
1195  } else {
1196  hobj->SetLineWidth(2);
1197  }
1198  if (oName.EndsWith("Mass")) hobj->Draw("e");
1199  else hobj->Draw();
1200  if (oName.BeginsWith("fcl_radial") && (hobj->GetEntries() > 0)) {
1201  //Fits to radial steps FTPCE+W/////05/14/08///nav+gvb
1202  hobj->GetXaxis()->SetRangeUser(6.5,9.0);
1203  hobj->Fit("pol0","","", 6.5, 7.2);
1204  double n1= hobj->GetFunction("pol0")->GetParameter(0);
1205  hobj->Fit("pol0","","",8.3,9.0);
1206  double n2= hobj->GetFunction("pol0")->GetParameter(0);
1207  fitFRS->SetParameters(n1, 7.85, 0.35, n2);
1208  hobj->Fit(fitFRS, "R");
1209 
1210  double rstep = fitFRS->GetParameter(1);
1211  double erstep = fitFRS->GetParError(1);
1212  float hmin = hobj->GetMinimum();
1213  float hmax = hobj->GetMaximum();
1214  ruler.SetLineColor(kBlack);
1215  ruler.SetLineWidth(2);
1216  ruler.DrawLine(7.8,hmin,7.8,hmax);
1217  ruler.SetLineColor(kGreen);
1218  ruler.SetLineWidth(3);
1219  ruler.DrawLine(rstep,hmin,rstep,hmax);
1220  ruler.SetLineWidth(1);
1221  ruler.DrawLine(rstep-erstep,hmin,rstep-erstep,hmax);
1222  ruler.DrawLine(rstep+erstep,hmin,rstep+erstep,hmax);
1223  }
1224 
1225 
1226  }
1227 
1228  if (oName.Contains("NullPrimVtx")) {
1229  int msdVtx = (int) (hobj->GetBinContent(hobj->FindBin(-1.)));
1230  int qstVtx = (int) (hobj->GetBinContent(hobj->FindBin(0.)));
1231  int goodVtx = (int) (hobj->GetBinContent(hobj->FindBin(1.)));
1232  int fndVtx = qstVtx + goodVtx;
1233  int totVtx = fndVtx + msdVtx;
1234  Float_t txtSiz = latex.GetTextSize();
1235  latex.SetTextSize(txtSiz*1.5);
1236  latex.SetTextAngle(90);
1237  latex.SetTextAlign(3);
1238  latex.SetTextColor(4);
1239  latex.DrawLatex(-0.8,0,Form(" missed: %d",msdVtx));
1240  latex.DrawLatex(0.2,0,Form(" questionable: %d",qstVtx));
1241  latex.DrawLatex(1.2,0,Form(" good: %d",goodVtx));
1242  latex.SetTextSize(txtSiz*2);
1243  latex.SetTextColor(2);
1244  latex.DrawLatex(-1.8,0,Form(" total: %d",totVtx));
1245  latex.SetTextAlign(1);
1246  latex.SetTextColor(kGreen+3);
1247  latex.DrawLatex(-1.1,0,Form(" found: %d",fndVtx));
1248  // restore
1249  latex.SetTextColor(1);
1250  latex.SetTextSize(txtSiz);
1251  LOG_INFO << (m_CurPrefix ? possiblePrefixes[m_CurPrefix] : "GE")
1252  << (analRepeat ? " Ref" : "")
1253  << " QA Events (found vtx/total) "
1254  << fndVtx << " / " << totVtx << endm;
1255  }
1256 
1257  if (oName.Contains("h1_inv_mass_cluster")) {
1258  ruler.SetLineColor(2);
1259  ruler.SetLineWidth(1);
1260  ruler.DrawLine(0.135,0.,0.135,hobj->GetMaximum());
1261  }
1262 
1263  if (oName.Contains("h2_cluster_dgg_vs_E1pE2")) {
1264  if (!fcsFunc1) {
1265  fcsFunc1= new TF1("Zgg=0","156.20/x",4,30);
1266  fcsFunc1->SetTitle("Z_{gg}=0");
1267  fcsFunc2= new TF1("Zgg=0.35","166.747/x",4,30);
1268  fcsFunc2->SetTitle("Z_{gg}=0.35");
1269  fcsFunc3= new TF1("Zgg=0.70","218.72/x",4,30);
1270  fcsFunc3->SetTitle("Z_{gg}=0.7");
1271  fcsFunc1->SetLineColor(2);
1272  fcsFunc2->SetLineColor(6);
1273  fcsFunc3->SetLineColor(28);
1274  }
1275  fcsFunc1->Draw("same");
1276  fcsFunc2->Draw("same");
1277  fcsFunc3->Draw("same");
1278  }
1279 
1280  if (oName.Contains("iTpcSector")) {
1281  // Draw RDO boundaries
1282  ruler.SetLineColor(1);
1283  ruler.SetLineWidth(1);
1284  // between RDOs 4-8, draw +/- (npads_row1+npads_row2)/2 * (pitch/2)
1285  float pitch = 0.67/2.0; // 6.7mm pitch
1286  ruler.DrawLine(-137*pitch,64.5,137*pitch,64.5);
1287  ruler.DrawLine(-123*pitch,56.5,123*pitch,56.5);
1288  ruler.DrawLine(-111*pitch,48.5,111*pitch,48.5);
1289  ruler.DrawLine( -97*pitch,40.5, 97*pitch,40.5);
1290  // between RDOs 1-2, outer X pads [X/2 at each end] are in RDO 1
1291  pitch = 0.50/2.0; // 5.0mm pitch
1292  int row_width = 70; int in_step1 = 52; int in_step2 = 54;
1293  float row1=10.5; float row2 = row1+1.0; float row3 = row1+2.0;
1294  ruler.DrawLine(-(row_width-in_step2)*pitch,row3,(row_width-in_step2)*pitch,row3);
1295  ruler.DrawLine(-(row_width-in_step2)*pitch,row2,-(row_width-in_step1)*pitch,row2);
1296  ruler.DrawLine((row_width-in_step2)*pitch,row2,(row_width-in_step1)*pitch,row2);
1297  ruler.DrawLine(-row_width*pitch,row1,-(row_width-in_step1)*pitch,row1);
1298  ruler.DrawLine((row_width-in_step1)*pitch,row1,row_width*pitch,row1);
1299  ruler.DrawLine(-(row_width-in_step1)*pitch,row1,-(row_width-in_step1)*pitch,row2);
1300  ruler.DrawLine((row_width-in_step1)*pitch,row1,(row_width-in_step1)*pitch,row2);
1301  ruler.DrawLine(-(row_width-in_step2)*pitch,row3,-(row_width-in_step2)*pitch,row2);
1302  ruler.DrawLine((row_width-in_step2)*pitch,row3,(row_width-in_step2)*pitch,row2);
1303  // between RDOs 2-3&4
1304  bool east = (atoi(&(oName.Data()[oName.Last('r')+1])) > 12);
1305  pitch = (east ? 0.50 : -0.50); // 5.0mm pitch, switch orientation for east vs. west
1306  float p01 = 46*pitch; float row01 = 22.5;
1307  float p02 = 36*pitch; float row02 = 22.5;
1308  float p03 = 36*pitch; float row03 = 21.5;
1309  float p04 = 39*pitch; float row04 = 21.5;
1310  float p05 = 39*pitch; float row05 = 20.5;
1311  float p06 = 19*pitch; float row06 = 20.5;
1312  float p07 = 19*pitch; float row07 = 22.5;
1313  float p08 = 9*pitch; float row08 = 22.5;
1314  float p09 = 9*pitch; float row09 = 23.5;
1315  float p10 = 2*pitch; float row10 = 23.5;
1316  float p11 = 2*pitch; float row11 = 24.5;
1317  ruler.DrawLine( p01,row01, p02,row02);
1318  ruler.DrawLine( p02,row02, p03,row03);
1319  ruler.DrawLine( p03,row03, p04,row04);
1320  ruler.DrawLine( p04,row04, p05,row05);
1321  ruler.DrawLine( p05,row05, p06,row06);
1322  ruler.DrawLine( p06,row06, p07,row07);
1323  ruler.DrawLine( p07,row07, p08,row08);
1324  ruler.DrawLine( p08,row08, p09,row09);
1325  ruler.DrawLine( p09,row09, p10,row10);
1326  ruler.DrawLine( p10,row10, p11,row11);
1327  ruler.DrawLine(-p01,row01,-p02,row02);
1328  ruler.DrawLine(-p02,row02,-p03,row03);
1329  ruler.DrawLine(-p03,row03,-p04,row04);
1330  ruler.DrawLine(-p04,row04,-p05,row05);
1331  ruler.DrawLine(-p05,row05,-p06,row06);
1332  ruler.DrawLine(-p06,row06,-p07,row07);
1333  ruler.DrawLine(-p07,row07,-p08,row08);
1334  ruler.DrawLine(-p08,row08,-p09,row09);
1335  ruler.DrawLine(-p09,row09,-p10,row10);
1336  ruler.DrawLine(-p10,row10,-p11,row11);
1337  ruler.DrawLine(-p11,row11, p11,row11);
1338  // between RDOs 3-4
1339  p01 = 0*pitch; row01 = 24.5;
1340  p02 = 0*pitch; row02 = 26.5;
1341  p03 = -7*pitch; row03 = 26.5;
1342  p04 = -7*pitch; row04 = 28.5;
1343  p05 = -8*pitch; row05 = 28.5;
1344  p06 = -8*pitch; row06 = 30.5;
1345  p07 = 10*pitch; row07 = 30.5;
1346  p08 = 10*pitch; row08 = 32.5;
1347  p09 = 11*pitch; row09 = 32.5;
1348  p10 = 11*pitch; row10 = 35.5;
1349  p11 = 10*pitch; row11 = 35.5;
1350  float p12 = 10*pitch; float row12 = 36.5;
1351  float p13 = 0*pitch; float row13 = 36.5;
1352  float p14 = 0*pitch; float row14 = 40.5;
1353  ruler.DrawLine(p01,row01,p02,row02);
1354  ruler.DrawLine(p02,row02,p03,row03);
1355  ruler.DrawLine(p03,row03,p04,row04);
1356  ruler.DrawLine(p04,row04,p05,row05);
1357  ruler.DrawLine(p05,row05,p06,row06);
1358  ruler.DrawLine(p06,row06,p07,row07);
1359  ruler.DrawLine(p07,row07,p08,row08);
1360  ruler.DrawLine(p08,row08,p09,row09);
1361  ruler.DrawLine(p09,row09,p10,row10);
1362  ruler.DrawLine(p10,row10,p11,row11);
1363  ruler.DrawLine(p11,row11,p12,row12);
1364  ruler.DrawLine(p12,row12,p13,row13);
1365  ruler.DrawLine(p13,row13,p14,row14);
1366  latex.SetTextAngle(0);
1367  latex.SetTextAlign(32);
1368  latex.DrawLatex(50,5,"RDO 1");
1369  latex.DrawLatex(50,18,"2");
1370  latex.DrawLatex(50,33,(east ? "4 , 3" : "3 , 4"));
1371  latex.DrawLatex(50,44,"5");
1372  latex.DrawLatex(50,52,"6");
1373  latex.DrawLatex(50,60,"7");
1374  latex.DrawLatex(50,68,"8");
1375 
1376  // Draw Anode guides
1377  ruler.SetLineColor(2);
1378  ruler.DrawLine(-52,10.5,-47,10.5);
1379  ruler.DrawLine(-52,20.5,-47,20.5);
1380  ruler.DrawLine(-52,30.5,-47,30.5);
1381  ruler.DrawLine(-52,40.5,-47,40.5);
1382  ruler.DrawLine(-52,48.1,-47,48.1);
1383  ruler.DrawLine(-52,56.1,-47,56.1);
1384  ruler.DrawLine(-52,64.1,-47,64.1);
1385  latex.SetTextAlign(12);
1386  latex.SetTextColor(2);
1387  latex.DrawLatex(-50, 5.1,"1 Anode");
1388  latex.DrawLatex(-50,15.1,"2");
1389  latex.DrawLatex(-50,25.1,"3");
1390  latex.DrawLatex(-50,35.1,"4");
1391  latex.DrawLatex(-50,44.1,"5");
1392  latex.DrawLatex(-50,52.1,"6");
1393  latex.DrawLatex(-50,60.1,"7");
1394  latex.DrawLatex(-50,68.1,"8");
1395  latex.SetTextColor(1);
1396  } else if (oName.Contains("TpcSector")) {
1397  // Draw RDO boundaries
1398  ruler.SetLineColor(1);
1399  ruler.SetLineWidth(1);
1400  // between RDOs 2-6, draw +/- (npads_row1+npads_row2)/2 * (pitch/2)
1401  float pitch = 0.67/2.0; // 6.7mm pitch
1402  ruler.DrawLine(-137*pitch,37.5,137*pitch,37.5);
1403  ruler.DrawLine(-123*pitch,29.5,123*pitch,29.5);
1404  ruler.DrawLine(-111*pitch,21.5,111*pitch,21.5);
1405  ruler.DrawLine( -97*pitch,13.5, 97*pitch,13.5);
1406  // between RDOs 1-2, outer 24(68) pads [12(34) at each end] are in RDO 1
1407  pitch = 0.335/2.0; // 3.35mm pitch
1408  int row_width = 134; int in_step = 68; float row1=6.5;
1409  //if (m_RunYear < 17) { row_width = 142; in_step = 24; row1=7.5; }
1410  float row2 = row1+1.0;
1411  ruler.DrawLine(-(row_width-in_step)*pitch,row2,(row_width-in_step)*pitch,row2);
1412  ruler.DrawLine(-row_width*pitch,row1,-(row_width-in_step)*pitch,row1);
1413  ruler.DrawLine((row_width-in_step)*pitch,row1,row_width*pitch,row1);
1414  ruler.DrawLine(-(row_width-in_step)*pitch,row1,-(row_width-in_step)*pitch,row2);
1415  ruler.DrawLine((row_width-in_step)*pitch,row1,(row_width-in_step)*pitch,row2);
1416  latex.SetTextAngle(0);
1417  latex.SetTextAlign(32);
1418  latex.DrawLatex(50,4,"RDO 1");
1419  latex.DrawLatex(50,10,"2");
1420  latex.DrawLatex(50,17,"3");
1421  latex.DrawLatex(50,25,"4");
1422  latex.DrawLatex(50,33,"5");
1423  latex.DrawLatex(50,41,"6");
1424 
1425  // Draw Anode guides
1426  ruler.SetLineColor(2);
1427  ruler.DrawLine(-52,2.9,-47,2.9);
1428  ruler.DrawLine(-52,6.2,-47,6.2);
1429  ruler.DrawLine(-52,9.4,-47,9.4);
1430  ruler.DrawLine(-52,13.5,-47,13.5);
1431  ruler.DrawLine(-52,21.1,-47,21.1);
1432  ruler.DrawLine(-52,29.1,-47,29.1);
1433  ruler.DrawLine(-52,37.1,-47,37.1);
1434  latex.SetTextAlign(12);
1435  latex.SetTextColor(2);
1436  latex.DrawLatex(-50,1.5,"1 Anode");
1437  latex.DrawLatex(-50,4.6,"2");
1438  latex.DrawLatex(-50,7.8,"3");
1439  latex.DrawLatex(-50,11.4,"4");
1440  latex.DrawLatex(-50,17.1,"5");
1441  latex.DrawLatex(-50,25.1,"6");
1442  latex.DrawLatex(-50,33.1,"7");
1443  latex.DrawLatex(-50,41.1,"8");
1444  latex.SetTextColor(1);
1445  }
1446 
1447  if (oName.Contains("PointRPTpc") ||
1448  oName.Contains("PointXYTpc") ) {
1449  // Draw sector boundaries and label
1450  int eastsec = (oName.Contains("TpcE") || oName.Contains("TpcQE") ? 12 : 0);
1451  ruler.SetLineColor(0);
1452  ruler.SetLineStyle(2);
1453  float sz = latex.GetTextSize();
1454  latex.SetTextAlign(32);
1455  latex.SetTextSize(0.032);
1456  int secn = (eastsec ? 21 : 13);
1457  float phistep = TMath::Pi()/12;
1458  for (float phi=phistep; phi<6.2; phi+=phistep) {
1459  float xsec = TMath::Cos(phi);
1460  float ysec = TMath::Sin(phi);
1461  ruler.DrawLine(50*xsec,50*ysec,199*xsec,199*ysec);
1462  phi+=phistep;
1463  xsec = TMath::Cos(phi);
1464  ysec = TMath::Sin(phi);
1465  latex.SetTextAngle(phi*180/TMath::Pi()+1); // +1 degree necessary to avoid PDF bug
1466  latex.DrawLatex(52*xsec,52*ysec,Form("%d",secn%12 + eastsec + 1));
1467  secn += (eastsec ? 1 : -1);
1468  }
1469  latex.SetTextSize(sz);
1470  ruler.SetLineStyle(1);
1471  }
1472 
1473  if (oName.EndsWith("QaPointPhiT")) {
1474  float hmin = hobj->GetMinimum() ;
1475  float hmax = hobj->GetMaximum() ;
1476  float hmid = hmin+0.12*(hmax-hmin);
1477  float sz = latex.GetTextSize();
1478  latex.SetTextSize(0.045);
1479  latex.SetTextColor(2);
1480  latex.SetTextAngle(90);
1481  latex.SetTextAlign(22);
1482  for (int secn = 1; secn < 13; secn++) {
1483  float phisec = 87 - secn*30;
1484  while (phisec < 1) phisec += 360;
1485  latex.DrawLatex(phisec,hmid,Form("%d / %d",secn,(23-secn)%12 + 13));
1486  }
1487  latex.SetTextSize(sz);
1488  latex.SetTextColor(1);
1489  latex.SetTextAngle(0);
1490  }
1491 
1492  if ((chkdim == 2) && (oName.Contains("TPC_charge"))) {
1493  TH2F* hobj_TPC_charge = (TH2F*) hobj;
1494 
1495  TProfile *prof_hobjTPC = (TProfile*) hobj_TPC_charge->ProfileX();
1496  prof_hobjTPC->SetMarkerStyle(20);
1497  prof_hobjTPC->SetMarkerSize(0.6);
1498  prof_hobjTPC->SetMarkerColor(6);
1499  prof_hobjTPC->SetLineColor(6);
1500 
1501  //landau fit
1502  static TF1 *flandau = 0;
1503  if (!flandau) flandau = new TF1("flandau", "[0]*TMath::Landau(x,[1],[2])",0,200);
1504  flandau->SetParameters(2e-3*(hobj->GetEntries()), 125, 30);
1505  hobj_TPC_charge->FitSlicesY(flandau);
1506  hobj->Draw("colz");
1507  TH1D *mean_landau = (TH1D*) (gDirectory->Get(Form("%s_1",hobj_TPC_charge->GetName())));
1508  if (mean_landau) {
1509  mean_landau->SetMarkerStyle(20);
1510  mean_landau->SetMarkerSize(0.6);
1511  mean_landau->SetMarkerColor(1);
1512  mean_landau->Draw("same p");
1513  } else {
1514  LOG_WARN << "TPC_charge Landau FitSlicesY failed for " << hobj->GetName() << endm;
1515  }
1516 
1517  //rms landau
1518  //TH1D *mean_landau2 = (TH1D*) (gDirectory->Get(Form("%s_2",hobj_TPC_charge_lan->GetName())));
1519  //if (mean_landau2) {
1520  // mean_landau2->SetMarkerStyle(20);
1521  // mean_landau2->SetMarkerSize(0.6);
1522  // mean_landau2->SetMarkerColor(2);
1523  // mean_landau2->Draw("same p");
1524  //}
1525 
1526  prof_hobjTPC->Draw("same p");
1527 
1528  //draw a legend
1529  TLegend *ilegend = new TLegend(0.15,0.80,0.38,0.90);
1530  ilegend->SetFillColor(0);
1531  ilegend->SetHeader("Legend");
1532  ilegend->SetMargin(0.25);
1533  ilegend->AddEntry(prof_hobjTPC,"Profile mean","pz");
1534  ilegend->AddEntry(mean_landau,"Landau MPV","pz");
1535  ilegend->Draw();
1536  }
1537 
1538  if (padAdvance) {if (gPad) gPad->Update();}
1539  else {if (!analRepeat) padCount--;}
1540 
1541  } // analRepeat for loop
1542 
1543  // Run reference analysis...(on first loop, before we forgot hobj)
1544  if (hobjR) {
1545  StHistUtilRef* huR = 0;
1546  if (m_refCuts) {
1547  // try: full name
1548  huR = (StHistUtilRef*) (m_refCuts->FindObject(oname));
1549  if (!huR) {
1550  // try: strip just maker from name
1551  Int_t tempint = -1;
1552  TString onamebase = StripPrefixes(oname,tempint,-1);
1553  huR = (StHistUtilRef*) (m_refCuts->FindObject(onamebase.Data()));
1554  if (!huR) {
1555  // try: strip just trigger type from name
1556  onamebase = StripPrefixes(oname,tempint,1);
1557  huR = (StHistUtilRef*) (m_refCuts->FindObject(onamebase.Data()));
1558  if (!huR) {
1559  // try: strip maker and trigger type from name
1560  onamebase = StripPrefixes(oname,tempint,0);
1561  huR = (StHistUtilRef*) (m_refCuts->FindObject(onamebase.Data()));
1562  }
1563  }
1564  }
1565  }
1566  double result = 0;
1567 
1568  // default to Kolm. max distance, no cut
1569  int mode = ( huR ? huR->Mode : 2 );
1570  double cut = ( huR ? huR->Cut : 0);
1571 
1572  // modes:
1573  switch (mode) {
1574  case (0) : // Chi2
1575  result = hobjO->Chi2Test(hobjR,(huR ? huR->Options() : "WW"));
1576  break;
1577  case (1) : // Kolmogorov probability
1578  // Note that Sumw2() seems to affect the way histograms are drawn,
1579  // but Kolmogorov test complains in current ROOT version (won't in
1580  // future versions). Don't know what to do for now.
1581  //if (hobjO->GetSumw2N() == 0) hobjO->Sumw2();
1582  //if (hobjR->GetSumw2N() == 0) hobjR->Sumw2();
1583  result = hobjO->KolmogorovTest(hobjR,(huR ? huR->Options() : ""));
1584  break;
1585  case (2) : // Kolmogorov maximum distance
1586  // Seems to work for both 1D and 2D
1587  TString analOpts = (huR ? huR->Options() : "");
1588  analOpts.ToUpper();
1589  if (! analOpts.Contains('M')) analOpts += 'M';
1590  result = 1.0 - hobjO->KolmogorovTest(hobjR,analOpts.Data());
1591  break;
1592  }
1593  bool analPass = (result >= cut);
1594 
1595  TString score = Form("Score: %4.2f (%s vs. %4.2f)",result,
1596  (analPass ? "PASS" : "FAIL"),cut);
1597 
1598  if (m_PrintMode & QAU1<<QAprintIndivRef)
1599  gPad->Print(Form("Ref_%s%s",oName.Data(),m_OutIndividuals.Data()));
1600  if (m_QAShiftMode)
1601  gPad->Print(Form("Ref_%s.png",oName.Data()));
1602 
1603  if (objPad) {
1604  objPad->cd();
1605  float sz = latex.GetTextSize();
1606  latex.SetTextSize(0.038);
1607  latex.SetTextAngle(90);
1608  latex.SetTextAlign(11);
1609  latex.SetTextColor(analPass ? kGreen+4 : kRed+2);
1610  latex.DrawTextNDC(0.999,0.001,score.Data());
1611  objPad->Update();
1612  latex.SetTextColor(1);
1613  latex.SetTextSize(sz);
1614  if (m_PrintMode & QAU1<<QAprintIndiv)
1615  gPad->Print(Form("%s%s",oName.Data(),m_OutIndividuals.Data()));
1616  if (m_QAShiftMode)
1617  gPad->Print(Form("%s.png",oName.Data()));
1618  }
1619 
1620 
1621  LOG_INFO << Form(" - %d. Reference Test (mode=%d) %s",
1622  histReadCounter,mode,score.Data()) << endm;
1623  if (strlen(m_refResultsFile)) {
1624  if (!R_ostr) R_ostr = new ofstream(m_refResultsFile);
1625  (*R_ostr) << m_CurPage << " " << curPad << " " << oName << " " << result << endl;
1626  }
1627  } else {
1628  if (m_PrintMode & QAU1<<QAprintIndiv)
1629  gPad->Print(Form("%s%s",oName.Data(),m_OutIndividuals.Data()));
1630  if (m_QAShiftMode)
1631  gPad->Print(Form("%s.png",oName.Data()));
1632  if (strlen(m_refResultsFile)) {
1633  if (!R_ostr) R_ostr = new ofstream(m_refResultsFile);
1634  (*R_ostr) << m_CurPage << " " << curPad << " " << oName << " 1.0" << endl;
1635  }
1636  }
1637 
1638  }
1639  }
1640 
1641 //NOTE! (13jan00,kt)
1642 //--> obj->Draw just appends the histogram to the list
1643 // --> you must update the current pad (gPad) or the whole big pad (graphPad)
1644 // to actually see the stupid thing
1645 
1646 // just ended actual loop over histograms !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1647  }
1648 
1649  if (C_ostr) delete C_ostr;
1650  if (R_ostr) delete R_ostr;
1651  if (root_ofile) delete root_ofile;
1652 
1653  CloseOutFile();
1654  return histCounter;
1655 }
1656 
1657 //_____________________________________________________________________________
1658 
1659 
1660 TList* StHistUtil::FindHists(const Char_t *dirName, const Char_t *withPrefix)
1661 {
1662 
1663 // NOTE - must have already used method SetPntrToMaker to get the
1664 // pointer m_PntrToMaker to an StMaker class!
1665 //
1666 
1667 // Method 1 ------------------------------------------------
1668 // Method FindHists -->
1669 // Find pointer to histograms under a Maker
1670 
1671  TList *dList=0;
1672 
1673  LOG_INFO << " Beg: FindHists, dList pointer = " << dList << endm;
1674 
1675 //---- First look under Maker for histograms ==>
1676 //They should show up in your Maker's directory, so search for them there,
1677 // i.e. MakerName/.hist is where they'd be
1678 // Note: Histograms is a method of StMaker
1679 //---- If you have a chain, you'll always have the .hist directory, so
1680 // have to check if there's really anything there (so use First method)
1681 
1682 //
1683  PathCopy(m_dirName,dirName);
1684  StMaker *temp = m_PntrToMaker->GetMaker(m_dirName);
1685 
1686  if (temp) {
1687  LOG_INFO << "FindHists - found pointer to maker" << endm;
1688  dList = temp->Histograms();
1689  }
1690 
1691 // Now check to see if any histograms exist here (look for something in
1692 // the list (test)
1693  TObject *test=0;
1694  if (dList) test = dList->First();
1695  if (test){
1696  LOG_INFO << " FindHists - found hist. in Maker-Branch " << endm;
1697  }
1698 
1699  LOG_INFO << " Mid: FindHists, dList pointer = " << dList << endm;
1700  LOG_INFO << " Mid: FindHists, test pointer = " << test << endm;
1701 
1702 // If you have the pointer but the hist. really aren't here, set
1703 // the pointer back to zero
1704  if (!test) dList = 0;
1705 
1706  LOG_INFO << " Mid2: FindHists, dList pointer = " << dList << endm;
1707  LOG_INFO << " Mid2: FindHists, test pointer = " << test << endm;
1708 
1709 
1710  if (!dList) {
1711 
1712 // Method 2 -----------------------------------------------------
1713 
1714 //-------------- Now try and see if they're in histBranch from output of bfc
1715 
1716 
1717  St_DataSet *hist=0;
1718  hist = m_PntrToMaker->GetDataSet("hist");
1719  if (hist) {
1720 // hist->ls(9);
1721 
1722 // must look in dirNameHist
1723 // use TString to append "Hist" to the dirName
1724 // += is overloaded operator of TString
1725 
1726  TString hBN(m_dirName);
1727  hBN += "Hist";
1728 
1729 //find particular branch
1730  St_DataSet *QAH = 0;
1731  QAH = hist->Find(hBN.Data());
1732 // or can create iterator and look over all branches
1733 
1734 //now get the list of histograms
1735  if (QAH) {
1736  dList = (TList *)QAH->GetObject();
1737  }
1738 
1739  }
1740 
1741 // now have we found them?
1742  if (dList){
1743  LOG_INFO << " FindHists - found hist. in histBranch, with name: "
1744  << m_dirName << endm;
1745  }
1746  else {
1747  LOG_INFO << " FindHists - histogram branch has not been found for branch --> "
1748  << m_dirName << endm;
1749  }
1750 
1751  }
1752 
1753  if (dList && (withPrefix || m_ListOfPrint)) dList = TrimListByPrefix(dList,withPrefix);
1754 
1755  LOG_INFO << " FindHists, dList pointer = " << dList << endm;
1756 
1757 
1758  return dList;
1759 }
1760 //_____________________________________________________________________________
1761 
1762 TList* StHistUtil::FindHists(TFile* histFile, const Char_t* withPrefix) {
1763  if (!histFile) return 0;
1764  TList* dList = histFile->GetList();
1765  if (dList->GetSize() == 0) {
1766  histFile->ReadAll();
1767  dList = histFile->GetList();
1768  }
1769 
1770  TObject* test = (dList ? dList->First() : 0);
1771  LOG_INFO << " Mid5: FindHists, dList pointer = " << dList << endm;
1772  LOG_INFO << " Mid5: FindHists, test pointer = " << test << endm;
1773  if (!test) dList = 0;
1774 
1775  if (dList && (withPrefix || m_ListOfPrint)) dList = TrimListByPrefix(dList,withPrefix);
1776  return dList;
1777 }
1778 //_____________________________________________________________________________
1779 
1780 TList* StHistUtil::TrimListByPrefix(TList* dList, const Char_t* withPrefix) {
1781  TList* dList2 = new TList;
1782 
1783  //Now want to loop over all histograms
1784  // Create an iterator
1785  TIter nextObj(dList);
1786  TObject *obj = 0;
1787  int withPrefixNumber = -1;
1788  int prefixNumber = -1;
1789  if (withPrefix) StripPrefixes(withPrefix,withPrefixNumber);
1790  while ((obj = nextObj())) {
1791  Bool_t addIt = kTRUE;
1792  if (withPrefix) {
1793  StripPrefixes(obj->GetName(),prefixNumber);
1794  if (prefixNumber != withPrefixNumber) addIt = kFALSE;
1795  }
1796  if (addIt && ((!m_ListOfPrint) ||
1797  (obj->InheritsFrom("TH1") &&
1798  m_ListOfPrint->FindObject(obj->GetName())))) dList2->AddLast(obj);
1799  }
1800  m_ItemsToClear->Add(dList2);
1801  return dList2;
1802 }
1803 //_____________________________________________________________________________
1804 
1805 
1806 Int_t StHistUtil::ListHists(const Char_t *dirName)
1807 {
1808 // Method ListHists -->
1809 // List of all histograms
1810 
1811  if (Debug()) {
1812  LOG_INFO << " **** Now in StHistUtil::ListHists **** " << endm;
1813  }
1814 
1815 // get the TList pointer to the histograms:
1816  TList* dirList = (m_PntrToMaker ? FindHists(dirName) : FindHists(m_PntrToPlainFile));
1817 
1818  if (!dirList) { LOG_INFO << " ListHists - histograms not available! " << endm; }
1819 
1820 //Now want to loop over all histograms
1821 // Create an iterator
1822  TIter nextObj(dirList);
1823  Int_t histReadCount = 0;
1824  TObject *obj = 0;
1825 
1826 // use = here instead of ==, because we are setting obj equal to nextObj and then seeing if it's T or F
1827  while ((obj = nextObj())) {
1828 
1829 // now check if obj is a histogram
1830  if (obj->InheritsFrom("TH1")) {
1831 
1832  histReadCount++;
1833 // \n means newline, \" means print a quote
1834 // LOG_INFO << Form(" %d. Have histogram Type %s, Name %s with Title=\"%s\"\n",histReadCount,obj->ClassName(),obj->GetName(),obj->GetTitle()) << endm;
1835  LOG_INFO << " ListHists: Hist No. " << histReadCount << ", Type: " << obj->ClassName()
1836  << ", Name: " << obj->GetName() << ", Title \"" << obj->GetTitle() << "\" "<< endm;
1837  }
1838  }
1839 
1840  LOG_INFO << " ListHists: Total No. Histograms Booked = " << histReadCount <<endm;
1841  return histReadCount;
1842 }
1843 
1844 
1845 //_____________________________________________________________________________
1846 
1847 Int_t StHistUtil::PrintInfoHists(TList *dirList, const Char_t *fname )
1848 {
1849 
1850  LOG_INFO << " **** Now in StHistUtil::PrintInfoHists **** " << endm;
1851  LOG_INFO << " output file = " << fname << endm;
1852 
1853  ofstream fout(fname);
1854 
1855  if (!dirList) { LOG_INFO << " PrintInfoHists - histograms not available! " << endm; }
1856 
1857  Int_t histInfoCount = 0;
1858 
1859  if (dirList){
1860 
1861 //Now want to loop over all histograms
1862 // Create an iterator
1863  TIter nextObj(dirList);
1864  TObject *obj = 0;
1865 
1866  LOG_INFO << " Hist #, Name, #Entries, Mean, RMS " << endm;
1867  fout << " Hist #, Name, #Entries, Mean, RMS " << endl;
1868 
1869 // use = instead of ==, because we are setting obj equal to nextObj and then seeing if it's T or F
1870 
1871  while ((obj = nextObj())) {
1872 
1873 // now check if obj is a histogram
1874  if (obj->InheritsFrom("TH1")) {
1875 
1876  histInfoCount++;
1877 
1878  LOG_INFO <<
1879  histInfoCount << " " <<
1880  obj->GetName() << " " <<
1881  ((TH1 *)obj)->GetEntries() << " " <<
1882  ((TH1 *)obj)->GetMean() << " " <<
1883  ((TH1 *)obj)->GetRMS() << " " <<
1884  endm;
1885 
1886  fout <<
1887  histInfoCount << " " <<
1888  obj->GetName() << " " <<
1889  ((TH1 *)obj)->GetEntries() << " " <<
1890  ((TH1 *)obj)->GetMean() << " " <<
1891  ((TH1 *)obj)->GetRMS() << " " <<
1892  endl;
1893 
1894  }
1895  }
1896  } // if dirList
1897 
1898  LOG_INFO << " PrintInfoHists: # hist read = " << histInfoCount <<endm;
1899 
1900  return histInfoCount;
1901 }
1902 
1903 
1904 //_____________________________________________________________________________
1905 
1906 
1907 Int_t StHistUtil::CopyHists(TList *dirList)
1908 {
1909  if (Debug()) {
1910  LOG_INFO << " **** Now in StHistUtil::CopyHists **** " << endm;
1911  }
1912 
1913  if (!dirList) { LOG_INFO << " StHistUtil::CopyHists - histogram Pointer not set! " << endm; }
1914 
1915 // create array of pointers to the new histograms I will create
1916 
1917  Int_t tempint,ijk=0;
1918  Int_t histCopyCount = 0;
1919 
1920  if (dirList){
1921  TIter nextObj(dirList);
1922  TObject *obj = 0;
1923  while ((obj = nextObj())) {
1924  if (obj->InheritsFrom("TH1") &&
1925  (!m_ListOfPrint || (m_ListOfPrint->FindObject(obj->GetName())))) {
1926  histCopyCount++;
1927  if (ijk >= maxHistCopy - 1){
1928  Int_t newMaxHistCopy = maxHistCopy * 4;
1929  TH1** temp1 = new TH1ptr[newMaxHistCopy];
1930  memset(temp1,0,newMaxHistCopy*sizeOfTH1Ptr);
1931  memcpy(temp1,newHist,maxHistCopy*sizeOfTH1Ptr);
1932  delete [] newHist;
1933  newHist = temp1;
1934  maxHistCopy = newMaxHistCopy;
1935  } // if ijk
1936  newHist[ijk] = ((TH1 *)obj->Clone());
1937  if (ignorePrefixes) {
1938  newHist[ijk]->SetName (StripPrefixes(newHist[ijk]->GetName (),tempint).Data());
1939  newHist[ijk]->SetTitle(StripPrefixes(newHist[ijk]->GetTitle(),tempint).Data());
1940  }
1941  ijk++;
1942  } // if obj
1943  } // while obj
1944  } // if dirList
1945 
1946  LOG_INFO << " ListHists: Total No. Histograms Copied = " <<
1947  histCopyCount <<endm;
1948 
1949 // Now see if we can find these copies:
1950  // Int_t imk = 0;
1951  //for (imk=0;imk<histCopyCount;imk++) {
1952  // if (newHist[imk]->InheritsFrom("TH1")) {
1953  // LOG_INFO << " !!! NEW Type: " << newHist[imk]->ClassName() <<
1954  // ", Name: " << newHist[imk]->GetName() <<
1955  // ", Title: " << newHist[imk]->GetTitle() <<
1956  // ", Max: " << ((TH1 *)newHist[imk])->GetMaximum() << endm;
1957  // }
1958  //}
1959 
1960  return histCopyCount;
1961 }
1962 
1963 //_____________________________________________________________________________
1964 
1965 // **** IMPORTANT! *****
1966 // THIS METHOD ASSUMES YOU HAVE ALREADY USED CopyHists TO PUT
1967 // HISTOGRAMS FROM 1 FILE INTO newHist array
1968 // -- this method is used in subsequent files!
1969 
1970 Int_t StHistUtil::AddHists(TList *dirList,Int_t numHistCopy)
1971 {
1972  if (Debug()) {
1973  LOG_INFO << " **** Now in StHistUtil::AddHists **** " << endm;
1974  }
1975  // LOG_INFO << " num hists to copy = " << numHistCopy << endm;
1976 
1977  if (!dirList) { LOG_INFO <<
1978  " StHistUtil::AddHists - histogram Pointer not set! " << endm; }
1979 
1980  Int_t histAddCount = 0;
1981  Int_t tempInt=0;
1982 
1983  if (dirList){
1984  if (numHistCopy < 0) numHistCopy = dirList->GetSize();
1985  TIter nextObj(dirList);
1986  TObject *obj = 0;
1987 
1988  while ((obj = nextObj())) {
1989  if (obj->InheritsFrom("TH1")) {
1990  TString oName = obj->GetName();
1991  if (ignorePrefixes) oName = StripPrefixes(oName.Data(),tempInt);
1992 // now want to add these histograms to the copied ones:
1993  Int_t tempint,imk = 0;
1994  Bool_t notfound = true;
1995  for (imk=0;imk<maxHistCopy;imk++) {
1996  if (newHist[imk]) {
1997  TString nName = newHist[imk]->GetName();
1998  if (ignorePrefixes) nName = StripPrefixes(nName.Data(),tempInt);
1999  if (! (nName.CompareTo(oName))) {
2000  //LOG_INFO << " ---- hist num to add --- " << imk << endm;
2001  newHist[imk]->Add((TH1 *)obj);
2002  histAddCount++;
2003  notfound = false;
2004  //LOG_INFO << " !!! Added histograms with Name: " << newHist[imk]->GetName() << endm;
2005  } // strcmp
2006  } else break; // ran out of existing hists
2007  } // loop over imk
2008 
2009  if (notfound) {
2010  // Hist doesn't exist yet, so let's add it
2011  if (imk >= maxHistCopy-1){
2012  Int_t newMaxHistCopy = maxHistCopy * 4;
2013  TH1** temp1 = new TH1ptr[newMaxHistCopy];
2014  memset(temp1,0,newMaxHistCopy*sizeOfTH1Ptr);
2015  memcpy(temp1,newHist,maxHistCopy*sizeOfTH1Ptr);
2016  delete [] newHist;
2017  newHist = temp1;
2018  maxHistCopy = newMaxHistCopy;
2019  } // if imk
2020  newHist[imk] = ((TH1 *)obj->Clone());
2021  if (ignorePrefixes) {
2022  newHist[imk]->SetName (StripPrefixes(newHist[imk]->GetName (),tempint).Data());
2023  newHist[imk]->SetTitle(StripPrefixes(newHist[imk]->GetTitle(),tempint).Data());
2024  }
2025  numHistCopy++;
2026  histAddCount++;
2027 
2028  } // notfound
2029  } // if obj inherits from th1
2030  } //while
2031  } //dirlist
2032 
2033  LOG_INFO << " StHistUtil::AddHists: Total No. Histograms Added = " <<
2034  histAddCount <<endm;
2035 
2036 
2037  return histAddCount;
2038 }
2039 
2040 //_____________________________________________________________________________
2041 
2042 
2043 Int_t StHistUtil::ExamineLogYList()
2044 {
2045 // Method ExamineLogYList
2046 // List of all histograms that will be drawn with logy scale
2047 
2048  if (Debug()) {
2049  LOG_INFO << " **** Now in StHistUtil::ExamineLogYList **** " << endm;
2050  }
2051 
2052 // m_ListOfLogY - is a list of log plots
2053 // construct a TObject
2054  TObject *obj = 0;
2055 // construct a TIter ==> () is an overloaded operator in TIter
2056  TIter nextObj(m_ListOfLogY);
2057  Int_t LogYCount = 0;
2058 
2059 // use = here instead of ==, because we are setting obj equal to nextObj and then seeing if it's T or F
2060  while ((obj = nextObj())) {
2061 
2062  if (Debug()) { LOG_INFO << " StHistUtil::ExamineLogYList has hist " << obj->GetName() << endm; }
2063  LogYCount++;
2064 
2065  }
2066 
2067  LOG_INFO << " Now in StHistUtil::ExamineLogYList, No. Hist. in LogY scale = " << LogYCount <<endm;
2068  return LogYCount;
2069 }
2070 
2071 //_____________________________________________________________________________
2072 
2073 
2074 Int_t StHistUtil::ExamineLogXList()
2075 {
2076 // Method ExamineLogXList
2077 // List of all histograms that will be drawn with logX scale
2078 
2079  if (Debug()) {
2080  LOG_INFO << " **** Now in StHistUtil::ExamineLogXList **** " << endm;
2081  }
2082 
2083 // m_ListOfLogX - is a list of log plots
2084 // construct a TObject
2085  TObject *obj = 0;
2086 // construct a TIter ==> () is an overloaded operator in TIter
2087  TIter nextObj(m_ListOfLogX);
2088  Int_t LogXCount = 0;
2089 
2090 // use = here instead of ==, because we are setting obj equal to nextObj and then seeing if it's T or F
2091  while ((obj = nextObj())) {
2092 
2093  if (Debug()) {
2094  LOG_INFO << " StHistUtil::ExamineLogXList has hist " << obj->GetName() << endm;
2095  }
2096  LogXCount++;
2097 
2098  }
2099 
2100  LOG_INFO << " Now in StHistUtil::ExamineLogXList, No. Hist. in LogX scale = " << LogXCount <<endm;
2101  return LogXCount;
2102 }
2103 
2104 //_____________________________________________________________________________
2105 
2106 
2107 Int_t StHistUtil::ExaminePrintList()
2108 {
2109 // Method ExaminePrintList
2110 // List of all histograms that will be drawn,printed
2111 
2112  if (Debug()) {
2113  LOG_INFO << " **** Now in StHistUtil::ExaminePrintList **** " << endm;
2114  }
2115 
2116 // m_ListOfPrint - is a list of hist to print,draw
2117 
2118 // check if there is a list
2119  if (!m_ListOfPrint){
2120  LOG_INFO << " no subset print list was setup - all hist in directory will be printed " << endm;
2121  // return PrintCount;
2122  return 0;
2123  }
2124 
2125 // construct a TObject
2126  TObject *obj = 0;
2127 
2128 // construct a TIter ==> () is an overloaded operator in TIter
2129  TIter nextObj(m_ListOfPrint);
2130  Int_t PrintCount = 0;
2131 
2132 // use = here instead of ==, because we are setting obj equal to nextObj and then seeing if it's T or F
2133  while ((obj = nextObj())) {
2134 
2135  if (Debug()) {
2136  LOG_INFO << " StHistUtil::ExaminePrintList has hist " << obj->GetName() << endm;
2137  }
2138  PrintCount++;
2139 
2140  }
2141 
2142  LOG_INFO << " Now in StHistUtil::ExaminePrintList, No. Hist. to Print,Draw = " << PrintCount <<endm;
2143  return m_ListOfPrint->GetSize();
2144 }
2145 
2146 //_____________________________________________________________________________
2147 
2148 
2149 Int_t StHistUtil::AddToLogYList(const Char_t *HistName){
2150 // Method AddToLogYList
2151 // making list of all histograms that we want drawn with LogY scale
2152 
2153  if (Debug()) {
2154  LOG_INFO << " **** Now in StHistUtil::AddToLogYList **** " << endm;
2155  }
2156 
2157 // Since I'm creating a new list, must delete it in the destructor!!
2158 //make a new TList on heap(persistant); have already defined m_ListOfLogY in header file
2159  if (!m_ListOfLogY) m_ListOfLogY = new TList;
2160 
2161 // the add method for TList requires a TObject input (also can use TObjString)
2162 // create TObjString on heap
2163  TObjString *HistNameObj = new TObjString(HistName);
2164 
2165 // - check if it's already on the list - use FindObject method of TList
2166  TObject *lobj = 0;
2167  lobj = m_ListOfLogY->FindObject(HistName);
2168 // now can use Add method of TList
2169  if (!lobj) {
2170  m_ListOfLogY->Add(HistNameObj);
2171  if (Debug()) {
2172  LOG_INFO << " StHistUtil::AddToLogYList: " << HistName <<endm;
2173  }
2174  } else {
2175  LOG_INFO << " StHistUtil::AddToLogYList: " << HistName << " already in list - not added" <<endm;
2176  }
2177 
2178 // return using a method of TList (inherits GetSize from TCollection)
2179  return m_ListOfLogY->GetSize();
2180 }
2181 
2182 
2183 //_____________________________________________________________________________
2184 
2185 
2186 Int_t StHistUtil::AddToLogXList(const Char_t *HistName){
2187 // Method AddToLogXList
2188 // making list of all histograms that we want drawn with LogX scale
2189 
2190  if (Debug()) {
2191  LOG_INFO << " **** Now in StHistUtil::AddToLogXList **** " << endm;
2192  }
2193 
2194 // Since I'm creating a new list, must delete it in the destructor!!
2195 //make a new TList on heap(persistant); have already defined m_ListOfLogX in header file
2196  if (!m_ListOfLogX) m_ListOfLogX = new TList;
2197 
2198 // the add method for TList requires a TObject input (also can use TObjString)
2199 // create TObjString on heap
2200  TObjString *HistNameObj = new TObjString(HistName);
2201 
2202 // - check if it's already on the list - use FindObject method of TList
2203  TObject *lobj = 0;
2204  lobj = m_ListOfLogX->FindObject(HistName);
2205 // now can use Add method of TList
2206  if (!lobj) {
2207  m_ListOfLogX->Add(HistNameObj);
2208  if (Debug()) {
2209  LOG_INFO << " StHistUtil::AddToLogXList: " << HistName <<endm;
2210  }
2211  } else {
2212  LOG_INFO << " StHistUtil::AddToLogXList: " << HistName << " already in list - not added" <<endm;
2213  }
2214 
2215 // return using a method of TList (inherits GetSize from TCollection)
2216  return m_ListOfLogX->GetSize();
2217 }
2218 
2219 
2220 //_____________________________________________________________________________
2221 
2222 
2223 Int_t StHistUtil::AddToPrintList(const Char_t *HistName){
2224 
2225 // Method AddToPrintList
2226 // making list of all histograms that we want drawn,printed
2227 
2228  if (Debug()) {
2229  LOG_INFO << " **** Now in StHistUtil::AddToPrintList **** " << endm;
2230  }
2231 
2232 // Since I'm creating a new list, must delete it in the destructor!!
2233 //make a new TList on heap(persistant); have already defined m_ListOfPrint in header file
2234  if (!m_ListOfPrint) m_ListOfPrint = new TList;
2235 
2236 // the add method for TList requires a TObject input (also can use TObjString)
2237 // create TObjString on heap
2238  TObjString *HistNameObj = new TObjString(HistName);
2239 
2240 // now can use Add method of TList
2241  if (!m_ListOfPrint->Contains(HistName)) {
2242  m_ListOfPrint->Add(HistNameObj);
2243  if (Debug()) {
2244  LOG_INFO << " StHistUtil::AddToPrintList: " << HistName <<endm;
2245  }
2246  } else {
2247  LOG_INFO << " StHistUtil::AddToPrintList: " << HistName << " already in list - not added" <<endm;
2248  }
2249 
2250 // return using a method of TList (inherits GetSize from TCollection)
2251  return m_ListOfPrint->GetSize();
2252 
2253 }
2254 
2255 //_____________________________________________________________________________
2256 
2257 Int_t StHistUtil::RemoveFromLogYList(const Char_t *HistName){
2258 // Method RemoveFromLogYList
2259 // remove hist from list that we want drawn with LogY scale
2260 
2261  if (Debug()) {
2262  LOG_INFO << " **** Now in StHistUtil::RemoveFromLogYList **** " << endm;
2263  }
2264 
2265 // check if list exists:
2266  if (m_ListOfLogY) {
2267 
2268 // the remove method for TList requires a TObject input
2269 // - check if it's on the list - use FindObject method of TList
2270  TObject *lobj = 0;
2271  lobj = m_ListOfLogY->FindObject(HistName);
2272 // now can use Remove method of TList
2273  if (lobj) {
2274  m_ListOfLogY->Remove(lobj);
2275  if (Debug()) {
2276  LOG_INFO << " RemoveLogYList: " << HistName << " has been removed from list" <<endm;
2277  }
2278  } else {
2279  LOG_INFO << " RemoveLogYList: " << HistName << " not on list - not removing" <<endm;
2280  }
2281 
2282  // return using a method of TList (inherits GetSize from TCollection)
2283  return m_ListOfLogY->GetSize();
2284  }
2285 
2286  LOG_INFO << " RemoveLogYList: " << HistName << " not on list - not removing" <<endm;
2287  return 0;
2288 }
2289 
2290 
2291 //_____________________________________________________________________________
2292 
2293 Int_t StHistUtil::RemoveFromLogXList(const Char_t *HistName){
2294 // Method RemoveFromLogXList
2295 // remove hist from list that we want drawn with LogX scale
2296 
2297  if (Debug()) {
2298  LOG_INFO << " **** Now in StHistUtil::RemoveFromLogXList **** " << endm;
2299  }
2300 
2301 // check if list exists:
2302  if (m_ListOfLogX) {
2303 
2304 // the remove method for TList requires a TObject input
2305 // - check if it's on the list - use FindObject method of TList
2306  TObject *lobj = 0;
2307  lobj = m_ListOfLogX->FindObject(HistName);
2308 // now can use Remove method of TList
2309  if (lobj) {
2310  m_ListOfLogX->Remove(lobj);
2311  if (Debug()) {
2312  LOG_INFO << " RemoveLogXList: " << HistName << " has been removed from list" <<endm;
2313  }
2314  } else {
2315  LOG_INFO << " RemoveLogXList: " << HistName << " not on list - not removing" <<endm;
2316  }
2317 
2318  // return using a method of TList (inherits GetSize from TCollection)
2319  return m_ListOfLogX->GetSize();
2320  }
2321 
2322  LOG_INFO << " RemoveLogXList: " << HistName << " not on list - not removing" <<endm;
2323  return 0;
2324 }
2325 
2326 
2327 //_____________________________________________________________________________
2328 
2329 Int_t StHistUtil::RemoveFromPrintList(const Char_t *HistName){
2330 // Method RemoveFromPrintList
2331 // remove hist from list that we want drawn,printed
2332 
2333  if (Debug()) {
2334  LOG_INFO << " **** Now in StHistUtil::RemoveFromPrintList **** " << endm;
2335  }
2336 
2337 // check if list exists:
2338  if (m_ListOfPrint) {
2339 
2340 // the remove method for TList requires a TObject input
2341 // - check if it's on the list - use FindObject method of TList
2342  TObject *lobj = 0;
2343  lobj = m_ListOfPrint->FindObject(HistName);
2344 // now can use Remove method of TList
2345  if (lobj) {
2346  m_ListOfPrint->Remove(lobj);
2347  if (Debug()) {
2348  LOG_INFO << " RemovePrintList: " << HistName << " has been removed from list" <<endm;
2349  }
2350  } else {
2351  LOG_INFO << " RemovePrintList: " << HistName << " not on list - not removing" <<endm;
2352  }
2353 
2354  // return using a method of TList (inherits GetSize from TCollection)
2355  return m_ListOfPrint->GetSize();
2356  }
2357 
2358  LOG_INFO << " RemovePrintList: " << HistName << " not on list - not removing" <<endm;
2359  return 0;
2360 }
2361 
2362 
2363 //_____________________________________________________________________________
2364 // Method SetDefaultLogYList
2365 // - create default list of histograms we want plotted in LogY scale
2366 
2367 void StHistUtil::SetDefaultLogYList(const Char_t *dirName)
2368 {
2369 // Method SetDefaultLogYList
2370 // - create default list of histograms we want plotted in LogY scale
2371 
2372  if (Debug()) {
2373  LOG_INFO << " **** Now in StHistUtil::SetDefaultLogYList **** " << endm;
2374  }
2375 
2376 
2377  PathCopy(m_dirName,dirName);
2378  TString type;
2379  if (!strcmp(m_dirName,"QA"))
2380  type = "Tab";
2381  else if (!strcmp(m_dirName,"EventQA"))
2382  type = "StE";
2383  else // default for now
2384  type = "StE";
2385 
2386  const Char_t* sdefList[] = {
2387  #include "St_QA_Maker/QAhlist_logy.h"
2388  };
2389 
2390  Int_t lengofList = sizeof(sdefList)/sizeOfCharPtr;
2391  Int_t numLog = 0;
2392  Int_t ilg = 0;
2393  for (ilg=0;ilg<lengofList;ilg++) {
2394  TString listString = sdefList[ilg];
2395  if (!(listString.BeginsWith("fcl") ||
2396  listString.BeginsWith("fms_qt_") ||
2397  listString.BeginsWith("fpd_channel_"))) {
2398  for (Int_t k=0; k<numOfPosPrefixes; k++) {
2399  ((listString = type) += possiblePrefixes[k]) += sdefList[ilg];
2400  numLog = AddToLogYList(listString.Data());
2401  }
2402  } else numLog = AddToLogYList(listString.Data());
2403  }
2404 
2405  LOG_INFO << " !!! StHistUtil::SetDefaultLogYList, # histogram put in list " << numLog << endm;
2406 
2407 }
2408 
2409 //_____________________________________________________________________________
2410 // Method SetDefaultLogXList
2411 // - create default list of histograms we want plotted in LogX scale
2412 
2413 void StHistUtil::SetDefaultLogXList(const Char_t *dirName)
2414 {
2415 // Method SetDefaultLogXList
2416 // - create default list of histograms we want plotted in LogX scale
2417 
2418  if (Debug()) {
2419  LOG_INFO << " **** Now in StHistUtil::SetDefaultLogXList **** " << endm;
2420  }
2421 
2422  PathCopy(m_dirName,dirName);
2423  TString type;
2424  if (!strcmp(m_dirName,"QA"))
2425  type = "Tab";
2426  else if (!strcmp(m_dirName,"EventQA"))
2427  type = "StE";
2428  else // default for now
2429  type = "StE";
2430 
2431  const Char_t* sdefList[] = {
2432  #include "St_QA_Maker/QAhlist_logx.h"
2433  };
2434 
2435  Int_t lengofList = sizeof(sdefList)/sizeOfCharPtr;
2436  Int_t numLog = 0;
2437  Int_t ilg = 0;
2438  for (ilg=0;ilg<lengofList;ilg++) {
2439  TString listString = sdefList[ilg];
2440  if (!(listString.BeginsWith("fcl") ||
2441  listString.BeginsWith("fms_qt_") ||
2442  listString.BeginsWith("fpd_channel_"))) {
2443  for (Int_t k=0; k<numOfPosPrefixes; k++) {
2444  ((listString = type) += possiblePrefixes[k]) += sdefList[ilg];
2445  numLog = AddToLogXList(listString.Data());
2446  }
2447  } else numLog = AddToLogXList(listString.Data());
2448  }
2449 
2450  LOG_INFO << " !!! StHistUtil::SetDefaultLogXList, # histogram put in list " << numLog << endm;
2451 
2452 }
2453 
2454 
2455 //_____________________________________________________________________________
2456 // Method SetDefaultPrintList
2457 // - create default list of histograms we want drawn,printed
2458 // - analType can be one of a number of predefined histogram lists,
2459 // or a file which has histogram names listed one per line
2460 // (files like $STAR/StRoot/St_QA_Maker/QAhlist*.h are properly
2461 // parsed as well, but histogram names are sufficient)
2462 
2463 void StHistUtil::SetDefaultPrintList(const Char_t *dirName, const Char_t *analType)
2464 {
2465 
2466  LOG_INFO << " **** Now in StHistUtil::SetDefaultPrintList **** " << endm;
2467 
2468  const Char_t **sdefList=0;
2469  Int_t lengofList = 0;
2470  bool mustDeleteList = false;
2471 
2472  PathCopy(m_dirName,dirName);
2473  TString type;
2474  if (!strcmp(m_dirName,"QA"))
2475  type = "Tab";
2476  else if (!strcmp(m_dirName,"EventQA"))
2477  type = "StE";
2478  else // default for now
2479  type = "StE";
2480 
2481 
2482 // If not analysis Type is set, then don't setup a list
2483  if ((!strcmp(analType,"")) || (!strcmp(analType,"All")) ) {
2484  LOG_INFO << " All histograms in directory will be printed/drawn, no list set" << endm;
2485  return;
2486  }
2487 
2488 // Cosmic Data Table QA list .................................................
2489  const Char_t* sdefList1[] = {
2490  #include "St_QA_Maker/QAhlist_QA_Cosmic.h"
2491  };
2492  if ((!strcmp(m_dirName,"QA")) && (!strcmp(analType,"Cosmic"))) {
2493  sdefList = sdefList1; lengofList = sizeof(sdefList1)/sizeOfCharPtr;
2494  }
2495 
2496 // Test Table QA list.........................................................
2497  const Char_t* sdefList2[] = {
2498  #include "St_QA_Maker/QAhlist_QA_TestQATable.h"
2499  };
2500  if ((!strcmp(m_dirName,"QA")) && (!strcmp(analType,"TestQATable"))) {
2501  sdefList = sdefList2; lengofList = sizeof(sdefList2)/sizeOfCharPtr;
2502  }
2503 
2504 // FTPC Table QA list.........................................................
2505  const Char_t* sdefList3[] = {
2506  #include "St_QA_Maker/QAhlist_QA_Ftpc.h"
2507  };
2508  if ((!strcmp(m_dirName,"QA")) && (!strcmp(analType,"Ftpc"))) {
2509  sdefList = sdefList3; lengofList = sizeof(sdefList3)/sizeOfCharPtr;
2510  }
2511 
2512 // FTPC Table QA list.........................................................
2513  const Char_t* sdefList4[] = {
2514  #include "St_QA_Maker/QAhlist_QA_MDC3.h"
2515  };
2516  if ((!strcmp(m_dirName,"FlowTag")) && (!strcmp(analType,"MDC3"))) {
2517  sdefList = sdefList4; lengofList = sizeof(sdefList4)/sizeOfCharPtr;
2518  }
2519 
2520 // St_QA_Maker histograms without svt and ftpc histograms.....................
2521  const Char_t* sdefList5[] = {
2522  #include "St_QA_Maker/QAhlist_QA_year1.h"
2523  };
2524  if ((!strcmp(m_dirName,"QA")) && (!strcmp(analType,"year1"))) {
2525  sdefList = sdefList5; lengofList = sizeof(sdefList5)/sizeOfCharPtr;
2526  }
2527 
2528 // St_QA_Maker histograms without the svt and ftpc histograms.................
2529  const Char_t* sdefList6[] = {
2530  #include "St_QA_Maker/QAhlist_EventQA_year1.h"
2531  };
2532  if ((!strcmp(m_dirName,"EventQA")) && (!strcmp(analType,"year1"))) {
2533  sdefList = sdefList6; lengofList = sizeof(sdefList6)/sizeOfCharPtr;
2534  }
2535 
2536 // St_QA_Maker histograms for QA shift........................................
2537  const Char_t* sdefList7[] = {
2538  #include "St_QA_Maker/QAhlist_QA_qa_shift.h"
2539  };
2540  if ((!strcmp(m_dirName,"QA")) && (!strcmp(analType,"qa_shift"))) {
2541  sdefList = sdefList7; lengofList = sizeof(sdefList7)/sizeOfCharPtr;
2542  }
2543 
2544 // St_QA_Maker histograms for QA shift........................................
2545  const Char_t* sdefList8[] = {
2546  #include "St_QA_Maker/QAhlist_EventQA_qa_shift.h"
2547  };
2548  if ((!strcmp(m_dirName,"EventQA")) && (!strcmp(analType,"qa_shift"))) {
2549  sdefList = sdefList8; lengofList = sizeof(sdefList8)/sizeOfCharPtr;
2550  }
2551 
2552 // St_QA_Maker histograms for tpcSectors......................................
2553  const Char_t* sdefList9[] = {
2554  #include "St_QA_Maker/QAhlist_tpcSectors.h"
2555  };
2556  if ((!strcmp(m_dirName,"EventQA")) && (!strcmp(analType,"tpcSectors"))) {
2557  sdefList = sdefList9; lengofList = sizeof(sdefList9)/sizeOfCharPtr;
2558  }
2559 
2560 // St_QA_Maker histograms for Svt.............................................
2561  const Char_t* sdefList10[] = {
2562  #include "St_QA_Maker/QAhlist_Svt.h"
2563  };
2564  if (!strcmp(analType,"Svt")) {
2565  sdefList = sdefList10; lengofList = sizeof(sdefList10)/sizeOfCharPtr;
2566  }
2567 
2568 // St_QA_Maker histograms for subsystems......................................
2569  const Char_t* sdefList11[] = {
2570  #include "St_QA_Maker/QAhlist_subsystems.h"
2571  };
2572  if (!strcmp(analType,"subsys")) {
2573  sdefList = sdefList11; lengofList = sizeof(sdefList11)/sizeOfCharPtr;
2574  }
2575 
2576  if (!sdefList) {
2577  // Try reading in a file as specified by analType
2578  ifstream analFile(analType);
2579  if (analFile.good()) {
2580  LOG_INFO << "Reading print list from: " << analType << endm;
2581  sdefList = new charptr[4096];
2582  mustDeleteList = true;
2583  char analBuffer[256];
2584  TString analString;
2585  Bool_t commenting = kFALSE;
2586  Int_t commentIdx = -1;
2587  while (analFile.getline(analBuffer,255)) {
2588  analString = analBuffer;
2589  if (commenting) {
2590  commentIdx = analString.Index("*/");
2591  if (commentIdx>=0) {
2592  analString.Remove(0,commentIdx+2);
2593  commenting = kFALSE;
2594  } else continue;
2595  }
2596  if (!commenting) {
2597  commentIdx = analString.Index("/*");
2598  if (commentIdx>=0) {
2599  analString.Remove(commentIdx);
2600  commenting = kTRUE;
2601  }
2602  }
2603  commentIdx = analString.Index("//");
2604  if (commentIdx>=0) analString.Remove(commentIdx);
2605  analString.Remove(TString::kBoth,' ').Remove(TString::kBoth,',').Remove(TString::kBoth,' ').Remove(TString::kBoth,'"');
2606  Int_t alen = analString.Length();
2607  if (alen<1) continue;
2608  Char_t* tempc = new Char_t[alen+1];
2609  strcpy(tempc,analString.Data());
2610  sdefList[lengofList] = tempc;
2611  lengofList++;
2612  }
2613  } else {
2614  LOG_WARN << "Could not find list of histograms specified by: " << analType << endm;
2615  }
2616  }
2617 
2618  Int_t numPrt = 0;
2619  Int_t ilg = 0;
2620  for (ilg=0;ilg<lengofList;ilg++) {
2621  TString ilgString = sdefList[ilg];
2622  if (mustDeleteList) delete [] sdefList[ilg];
2623  Bool_t addIt = kTRUE;
2624  if (ilgString.BeginsWith(":")) {
2625  Ssiz_t endDetSpec = ilgString.Index(":",1) + 1;
2626  TString detSpec = ilgString(0,endDetSpec);
2627  addIt = DetectorIn(detSpec.Data());
2628  if (addIt) ilgString.Remove(0,endDetSpec);
2629  }
2630  if (addIt) {
2631  TRegexp fcsPattern("^h[12]_.+"); Ssiz_t ext=0;
2632  if (!(ilgString.BeginsWith("fcl") ||
2633  ilgString.BeginsWith("fms_qt_") ||
2634  ilgString.BeginsWith("Z3A") ||
2635  ilgString.Contains("etof",TString::ECaseCompare::kIgnoreCase) ||
2636  ilgString.Contains("_matchCand_") ||
2637  ilgString.BeginsWith("G_primary") ||
2638  ilgString.Index(fcsPattern,&ext)==0 ||
2639  ilgString.BeginsWith("fpd_channel_"))) {
2640  for (Int_t k=0; k<numOfPosPrefixes; k++) {
2641  TString listString = type;
2642  (listString += possiblePrefixes[k]) += ilgString;
2643  numPrt = AddToPrintList(listString.Data());
2644  }
2645  } else numPrt = AddToPrintList(ilgString.Data());
2646  if (Debug()) {
2647  LOG_INFO << " !!! adding histogram " << ilgString << " to print list " << endm ;
2648  }
2649  }
2650  }
2651  if (mustDeleteList) delete [] sdefList;
2652 
2653  LOG_INFO << " !!! StHistUtil::SetDefaultPrintList, # histogram put in list " << numPrt << endm;
2654 
2655 }
2656 
2657 //_____________________________________________________________________________
2658 
2659 // Method Overlay1D
2660 // - takes two TH1F histograms and overlays them
2661 
2662 Int_t StHistUtil::Overlay1D(Char_t *dirName,Char_t *inHist1,
2663  Char_t *inHist2) {
2664 
2665  LOG_INFO << " **** Now in StHistUtil::Overlay1D **** " << endm;
2666 
2667  Int_t n1dHists = 0;
2668  PathCopy(m_dirName,dirName);
2669 
2670 // get the TList pointer to the histograms
2671  TList* dirList = (m_PntrToMaker ? FindHists(m_dirName) : FindHists(m_PntrToPlainFile));
2672 
2673 // check that directory exists
2674  if (!dirList)
2675  return kStErr;
2676 
2677  LOG_INFO << "Histogram directory exists -> Find and overlay histograms" << endm;
2678 // Now want to loop over all histograms
2679 // Create an iterator
2680  TIter nextObj(dirList);
2681  TObject *obj = 0;
2682 
2683 // temporary holder histograms
2684  TH1F *hist1f1 = new TH1F;
2685  TH1F *hist1f2 = new TH1F;
2686 
2687 // use = here instead of ==, because we are setting obj equal to nextObj
2688 // and then seeing if it's T or F
2689  while ((obj = nextObj())) {
2690 
2691 // now check if obj is a histogram and see if it matches input name
2692  if (obj->InheritsFrom("TH1")) {
2693  if (obj->GetName() == (TString)inHist1 ||
2694  obj->GetName() == (TString)inHist2) {
2695  LOG_INFO << " Found Histogram: Type '" << obj->ClassName() << "', Name '"
2696  << obj->GetName() << "', Title '" << obj->GetTitle() << "'"
2697  << endm;
2698 
2699 // check on type of histogram and make copies
2700  if (obj->ClassName() == (TString)"TH1F") {
2701  if (obj->GetName() == (TString)inHist1) {
2702  *hist1f1 = *(TH1F *)obj;
2703  n1dHists++;
2704  }
2705  if (obj->GetName() == (TString)inHist2) {
2706  *hist1f2 = *(TH1F *)obj;
2707  n1dHists++;
2708  }
2709  } else {
2710  LOG_INFO << " ERROR: histogram not of type TH1F !!!" << endm;
2711  }
2712  }
2713  }
2714  }
2715 
2716 // if the two histograms exist, overlay them
2717  if (n1dHists == 2) {
2718  hist1f1->SetLineColor(4);
2719  hist1f1->SetLineStyle(1);
2720  hist1f2->SetLineColor(2);
2721  hist1f2->SetLineStyle(2);
2722 
2723  hist1f1->SetStats(kFALSE);
2724  hist1f2->SetStats(kFALSE);
2725 
2726  hist1f1->SetTitle(hist1f1->GetTitle()+(TString)" and "+hist1f2->GetTitle());
2727  hist1f2->SetTitle(hist1f1->GetTitle());
2728 // create a new canvas
2729  TCanvas *newCanvas = new TCanvas("c1d","Combined 1D Histogram",600,780);
2730  newCanvas->Draw();
2731 
2732 // write title at top of canvas
2733  Ltitle = new TPaveLabel(0.1,0.96,0.9,1.0,m_GlobalTitle.Data(),"br");
2734  Ltitle->SetFillColor(18);
2735  Ltitle->SetTextFont(32);
2736  Ltitle->SetTextSize(0.5);
2737  Ltitle->Draw();
2738 
2739 // now put in date & time at bottom right of canvas
2740  TDatime HistTime;
2741  const Char_t *myTime = HistTime.AsString();
2742  TPaveLabel *Ldatetime = new TPaveLabel(0.7,0.01,0.95,0.03,myTime,"br");
2743  Ldatetime->SetTextSize(0.6);
2744  Ldatetime->Draw();
2745 
2746 // create a pad
2747  TPad *newPad = new TPad("p1d","Combined 1D Histogram",0.02,0.04,0.98,0.93);
2748  newPad->Draw();
2749  newPad->cd();
2750 
2751 // draw the histograms
2752  if (hist1f1->GetMaximum() >= hist1f2->GetMaximum()) {
2753  hist1f1->Draw();
2754  hist1f2->Draw("Same");
2755  }
2756  else {
2757  hist1f2->Draw();
2758  hist1f1->Draw("Same");
2759  }
2760 
2761 // make a legend
2762  TLegend *legend = new TLegend(0.75,0.85,0.98,0.95);
2763  legend->SetFillColor(0);
2764  legend->SetHeader("Legend");
2765  legend->SetMargin(0.25);
2766  legend->AddEntry(hist1f1,inHist1,"l");
2767  legend->AddEntry(hist1f2,inHist2,"l");
2768  legend->Draw();
2769 
2770  newCanvas->Update();
2771 
2772  return kStOk;
2773  }
2774 
2775  return kStErr;
2776 }
2777 
2778 //_____________________________________________________________________________
2779 
2780 // Method Overlay2D
2781 // - takes two TH2F histograms and overlays them
2782 
2783 Int_t StHistUtil::Overlay2D(Char_t *dirName,Char_t *inHist1,
2784  Char_t *inHist2) {
2785 
2786  LOG_INFO << " **** Now in StHistUtil::Overlay2D **** " << endm;
2787 
2788  Int_t n2dHists = 0;
2789  PathCopy(m_dirName,dirName);
2790 
2791 // get the TList pointer to the histograms
2792  TList* dirList = (m_PntrToMaker ? FindHists(m_dirName) : FindHists(m_PntrToPlainFile));
2793 
2794 // check that directory exists
2795  if (!dirList)
2796  return kStErr;
2797 
2798  LOG_INFO << "Histogram directory exists -> Find and overlay histograms" << endm;
2799 // Now want to loop over all histograms
2800 // Create an iterator
2801  TIter nextObj(dirList);
2802  TObject *obj = 0;
2803 
2804 // temporary holder histograms
2805  TH2F *hist2f1 = new TH2F;
2806  TH2F *hist2f2 = new TH2F;
2807 
2808 // use = here instead of ==, because we are setting obj equal to nextObj
2809 // and then seeing if it's T or F
2810  while ((obj = nextObj())) {
2811 
2812 // now check if obj is a histogram and see if it matches input name
2813  if (obj->InheritsFrom("TH1")) {
2814  if (obj->GetName() == (TString)inHist1 ||
2815  obj->GetName() == (TString)inHist2) {
2816  LOG_INFO << " Found Histogram: Type '" << obj->ClassName() << "', Name '"
2817  << obj->GetName() << "', Title '" << obj->GetTitle() << "'"
2818  << endm;
2819 
2820 // check on type of histogram and make copies
2821  if (obj->ClassName() == (TString)"TH2F") {
2822  if (obj->GetName() == (TString)inHist1) {
2823  *hist2f1 = *(TH2F *)obj;
2824  n2dHists++;
2825  }
2826  if (obj->GetName() == (TString)inHist2) {
2827  *hist2f2 = *(TH2F *)obj;
2828  n2dHists++;
2829  }
2830  } else {
2831  LOG_INFO << " ERROR: histogram is not of type TH2F !!!" << endm;
2832  }
2833  }
2834  }
2835  }
2836 
2837 // if the two histograms exist, overlay them
2838  if (n2dHists == 2) {
2839  hist2f1->SetLineColor(4);
2840  hist2f2->SetLineColor(2);
2841 
2842  hist2f1->SetStats(kFALSE);
2843  hist2f2->SetStats(kFALSE);
2844 
2845  hist2f1->SetTitle(hist2f1->GetTitle()+(TString)" and "+hist2f2->GetTitle());
2846  hist2f2->SetTitle(hist2f1->GetTitle());
2847 
2848 // create a new canvas and pad to write to
2849  TCanvas *newCanvas = new TCanvas("c2d","Combined 2D Histogram",600,780);
2850  newCanvas->Draw();
2851 
2852 // write title at top of canvas
2853  Ltitle = new TPaveLabel(0.1,0.96,0.9,1.0,m_GlobalTitle.Data(),"br");
2854  Ltitle->SetFillColor(18);
2855  Ltitle->SetTextFont(32);
2856  Ltitle->SetTextSize(0.5);
2857  Ltitle->Draw();
2858 
2859 // now put in date & time at bottom right of canvas
2860  TDatime HistTime;
2861  const Char_t *myTime = HistTime.AsString();
2862  TPaveLabel *Ldatetime = new TPaveLabel(0.7,0.01,0.95,0.03,myTime,"br");
2863  Ldatetime->SetTextSize(0.6);
2864  Ldatetime->Draw();
2865 
2866 // create a pad
2867  TPad *newPad = new TPad("p2d","Combined 2D Histogram",0.02,0.04,0.98,0.93);
2868  newPad->Draw();
2869  newPad->cd();
2870 
2871 // draw the histograms
2872  if (hist2f1->GetMaximum() >= hist2f2->GetMaximum()) {
2873  hist2f1->Draw("Box");
2874  hist2f2->Draw("BoxSame");
2875  }
2876  else {
2877  hist2f2->Draw("Box");
2878  hist2f1->Draw("BoxSame");
2879  }
2880 
2881 // make a legend
2882  TLegend *legend = new TLegend(0.75,0.85,0.98,0.95);
2883  legend->SetFillColor(0);
2884  legend->SetHeader("Legend");
2885  legend->SetMargin(0.25);
2886  legend->AddEntry(hist2f1,inHist1,"f");
2887  legend->AddEntry(hist2f2,inHist2,"f");
2888  legend->Draw();
2889 
2890  newCanvas->Update();
2891 
2892  return kStOk;
2893  }
2894 
2895  return kStErr;
2896 }
2897 
2898 //_____________________________________________________________________________
2899 
2900 // Method GetRunYear
2901 // - determines the run year from the filename
2902 // - assumes runnumber is first all digit token after 2 "_" delimiters
2903 // - assumes runyear is all but the last 6 digits of runnumber
2904 
2905 Int_t StHistUtil::GetRunYear(const Char_t *filename) {
2906 
2907  m_RunYear = 0;
2908  TString FileName = filename;
2909  TString delim = "_";
2910  TObjArray* tokens = FileName.Tokenize(delim);
2911  for (int tk=2; tk<tokens->GetEntries(); tk++) {
2912  TString& tok = ((TObjString*) (tokens->At(tk)))->String();
2913  if (tok.IsDigit()) {
2914  Ssiz_t loc = tok.Length()-6;
2915  if (loc>0) m_RunYear = atoi(tok.Remove(loc).Data());
2916  break;
2917  }
2918  }
2919  return m_RunYear;
2920 
2921 }
2922 
2923 //_____________________________________________________________________________
2924 
2925 // Method SetRefAnalysis
2926 // - Properly set input and output files for reference histogram analysis
2927 
2928 void StHistUtil::SetRefAnalysis(const Char_t* refOutFile, const Char_t* refResultsFile,
2929  const Char_t* refCutsFile, const Char_t* refInFile) {
2930 
2931  LOG_INFO << "StHistUtil: Will run in reference histogram analysis mode." << endm;
2932  m_analMode = kTRUE;
2933 
2934  if (refInFile && strlen(refInFile)) {
2935  LOG_INFO << "StHistUtil: Using reference histogram file " << refInFile << endm;
2936  m_refInFile = new TFile(refInFile);
2937  if (!m_refInFile) {LOG_ERROR << "file not found: " << refInFile << endm;}
2938  }
2939  if (refCutsFile && strlen(refCutsFile)) {
2940  LOG_INFO << "StHistUtil: Using reference cuts file " << refCutsFile << endm;
2941  m_refCuts = new TList;
2942  ifstream refCuts(refCutsFile);
2943  if (!refCuts.is_open()) {
2944  LOG_ERROR << "StHistUtil: Unable to open cuts file! Proceeding with no cuts..." << endm;
2945  return;
2946  }
2947  char buf_name[256];
2948  char buf_opts[64];
2949  while (refCuts.good()) {
2950  refCuts >> buf_name;
2951  if (!refCuts.good()) break;
2952  int mode;
2953  double cut;
2954  refCuts >> mode >> cut >> buf_opts;
2955  if (!refCuts.good()) break;
2956  if (buf_opts[0] == '!') buf_opts[0] = 0; // no options
2957  if (Debug()) {
2958  LOG_INFO << "StHistUtil: Loading cut : " << buf_name << " : " << mode
2959  << " : " << cut << " : " << buf_opts << " :" << endm;
2960  }
2961  m_refCuts->AddLast(new StHistUtilRef(buf_name,buf_opts,mode,cut));
2962  }
2963  }
2964 
2965  // refOutFile will not be used if no reference histograms are found
2966  PathCopy(m_refResultsFile,refResultsFile);
2967  // refOutFile will not be used if already writing hists to a ROOT file
2968  PathCopy(m_refOutFile,refOutFile);
2969 }
2970 
2971 //_____________________________________________________________________________
2972 
2973 void StHistUtil::SetDetectors(const Char_t *detectors) {
2974  // entering a null pointer uses all detector hists
2975  // entering a zero length string uses no detector hists
2976  if (detectors) {
2977  m_Detectors = detectors;
2978  // use colons as delimiters
2979  m_Detectors.ReplaceAll(" ",":");
2980  m_Detectors.ReplaceAll(",",":");
2981  m_Detectors.Prepend(":");
2982  m_Detectors.Append(":");
2983  while (m_Detectors.Index("::") >= 0) m_Detectors.ReplaceAll("::",":");
2984  LOG_INFO << "StHistUtil::SetDetectors(): using detectors " << m_Detectors << endm;
2985  } else {
2986  m_Detectors = "";
2987  LOG_INFO << "StHistUtil::SetDetectors(): using all detectors" << endm;
2988  }
2989 }
2990 
2991 //_____________________________________________________________________________
2992 
2993 Bool_t StHistUtil::DetectorIn(const Char_t *detector) {
2994  // detector should be passed with delimeters, e.g. ":tpc:"
2995  // semicolons should separate multiple possible detectors ("or")
2996  // commas should separate multiple required detectors ("and")
2997  // Example ":tpc,svt;tpx,svt:" is either tpc+svt or tpx+svt
2998  Bool_t isIn = kFALSE;
2999  if (m_Detectors.Length()) {
3000  Ssiz_t idx = 0;
3001  TString detOpts = detector;
3002  TString curDetSet = "";
3003  TString curDet = "";
3004  while (idx < detOpts.Length() && ! isIn) {
3005  Ssiz_t D1 = detOpts.Index(";",idx);
3006  Ssiz_t setLength = (D1 < idx ? detOpts.Length() : D1+1) - idx;
3007  curDetSet = detOpts(idx,setLength);
3008  Ssiz_t idx2 = 0;
3009  isIn = kTRUE;
3010  while (idx2 < curDetSet.Length() && isIn) {
3011  Ssiz_t D2 = curDetSet.Index(",",idx2);
3012  Ssiz_t detLength = (D2 < idx2 ? curDetSet.Length() : D2+1) - idx2;
3013  curDet = curDetSet(idx2,detLength);
3014  curDet.ReplaceAll(";",":");
3015  curDet.ReplaceAll(",",":");
3016  curDet.Prepend(":");
3017  curDet.Append(":");
3018  while (curDet.Index("::") >= 0) curDet.ReplaceAll("::",":");
3019  if (m_Detectors.Index(curDet.Data()) < 0) isIn = kFALSE;
3020  else idx2 += detLength;
3021  }
3022  idx += setLength;
3023  }
3024  } else {
3025  // set to use all detectors
3026  isIn = kTRUE;
3027  }
3028  return isIn;
3029 }
3030 
3031 //_____________________________________________________________________________
3032 
3033 TH1* StHistUtil::FlipAxes(TH1* hist) {
3034  if (!(hist->InheritsFrom("TH2"))) return hist;
3035  // For now, just converts whatever into TH2D
3036  // ...and there is no automatic garbage collection
3037  Int_t xbins = hist->GetNbinsX();
3038  Int_t ybins = hist->GetNbinsY();
3039  TH2D* newhist = new TH2D(Form("%s_flip",hist->GetName()),hist->GetTitle(),
3040  ybins,hist->GetYaxis()->GetXmin(),hist->GetYaxis()->GetXmax(),
3041  xbins,hist->GetXaxis()->GetXmin(),hist->GetXaxis()->GetXmax());
3042  for (int xbin=1; xbin<=xbins; xbin++) {
3043  for (int ybin=1; ybin<=ybins; ybin++) {
3044  newhist->SetBinContent(ybin,xbin,hist->GetBinContent(xbin,ybin));
3045  newhist->SetBinError(ybin,xbin,hist->GetBinError(xbin,ybin));
3046  }
3047  }
3048  newhist->SetMinimum(hist->GetMinimumStored());
3049  newhist->SetMaximum(hist->GetMaximumStored());
3050  return newhist;
3051 }
3052 
3053 //_____________________________________________________________________________
3054 
3055 void StHistUtil::PathCopy(char* destination, const char* source) {
3056  // carefully copy path strings to avoid buffer over-writes or over-reads
3057  if (source && strcmp(destination,source)) {
3058  if (strlen(source) < maxPathLen) {
3059  strncpy(destination,source,maxPathLen-1);
3060  destination[maxPathLen-1] = 0;
3061  } else {
3062  LOG_ERROR << " source path too long: " << source << endm;
3063  }
3064  }
3065 }
3066 
3067 //_____________________________________________________________________________
3068 
3069 //_____________________________________________________________________________
3070 
3071 
3072 ClassImp(StHistUtilRef)
3073 StHistUtilRef::StHistUtilRef(const char* name, const char* opts,
3074  const int mode, const double cut):
3075  TNamed(name,opts),Mode(mode),Cut(cut) {}
Definition: IO.h:20
Definition: Cut.h:18
Definition: QAH.h:34
virtual TObject * GetObject() const
The depricated method (left here for the sake of the backward compatibility)
Definition: TDataSet.cxx:428
Definition: Stypes.h:44
Definition: Stypes.h:41
virtual TDataSet * Find(const char *path) const
Definition: TDataSet.cxx:362