Jpp - the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JGizmoToolkit.hh
Go to the documentation of this file.
1 #ifndef __JGIZMOTOOLKIT__
2 #define __JGIZMOTOOLKIT__
3 
4 #include <string>
5 #include <map>
6 #include <cmath>
7 
8 #include "TError.h"
9 #include "TFile.h"
10 #include "TClass.h"
11 #include "TObject.h"
12 #include "TKey.h"
13 #include "TH1.h"
14 #include "TH2.h"
15 #include "TGraph.h"
16 #include "TGraph2D.h"
17 #include "TMultiGraph.h"
18 #include "TString.h"
19 #include "TRegexp.h"
20 #include "TFormula.h"
21 #include "TF1.h"
22 #include "TF2.h"
23 #include "TIterator.h"
24 #include "TMethod.h"
25 #include "TMethodCall.h"
26 #include "TAxis.h"
27 #include "TMath.h"
28 
29 #include "JLang/JException.hh"
30 #include "JGizmo/JRootObjectID.hh"
31 
33 
34 
35 /**
36  * \author mdejong
37  */
38 
39 namespace JGIZMO {}
40 namespace JPP { using namespace JGIZMO; }
41 
42 /**
43  * Auxiliary applications for use of ROOT and more.
44  */
45 namespace JGIZMO {
46 
47  using JLANG::JParseError;
50 
51 
52  /**
53  * Auxiliary data structure for JOpera1D.cc and JOpera2D.cc applications.
54  */
55  struct JOpera {
56  //
57  // Histogram name.
58  //
59  static const char* const SAME_AS_OPERATION() { return "%"; } //!< Set name of output histogram to name of operation
60  static const char* const SAME_AS_INPUT() { return "="; } //!< Set name of output histogram to name of input histogram
61 
62  //
63  // Histogram operations.
64  //
65  static const char* const Add() { return "Add"; } //!< ROOT TH1::Add
66  static const char* const add() { return "add"; } //!< Add contents with lookup bin in second histogram
67  static const char* const Subtract() { return "Subtract"; } //!< ROOT TH1::Subtract
68  static const char* const subtract() { return "subtract"; } //!< Subtract contents with lookup bin in second histogram
69  static const char* const Multiply() { return "Multiply"; } //!< ROOT TH1::Multiply
70  static const char* const multiply() { return "multiply"; } //!< Multiply contents with lookup bin in second histogram
71  static const char* const Divide() { return "Divide"; } //!< ROOT TH1::Divide
72  static const char* const divide() { return "divide"; } //!< Divide contents with lookup bin in second histogram
73  static const char* const efficiency() { return "efficiency"; } //!< Divide contents and multiply errors with inefficiency
74  static const char* const stdev() { return "stdev"; } //!< Set contents to standard deviation
75  static const char* const sqrt() { return "sqrt"; } //!< Set contents to signed difference between squares
76  static const char* const Replace() { return "Replace"; } //!< Set contents to associated function
77  };
78 
79 
80  /**
81  * Time stamp of earliest UTC time.
82  */
83  static const char* const TIMESTAMP = "#splitline{}{#splitline{%d:%m:%y}{ %H:%M}}%F1970-01-01 00:00:00";
84 
85 
86  /**
87  * Get TFile pointer corresponding to give file name.
88  *
89  * The TFile pointer of an already opened file is recovered, else a new file is opened.\n
90  * Note that the closure of the opened files should be done by the caller of this method.
91  *
92  * \param file_name file name
93  * \param option TFile::Open option
94  * \return pointer to TFile
95  */
96  inline TFile* getFile(const std::string& file_name, const std::string& option = "exist")
97  {
98  using namespace std;
99 
100  gErrorIgnoreLevel = kError;
101 
102  static map<string, TFile*> zmap;
103 
104  map<string, TFile*>::iterator i = zmap.find(file_name);
105 
106  if (i == zmap.end() || i->second == NULL || !i->second->IsOpen()) {
107 
108  TFile* file = TFile::Open(file_name.c_str(), option.c_str());
109 
110  zmap[file_name] = file;
111 
112  return file;
113 
114  } else {
115 
116  return i->second;
117  }
118  }
119 
120 
121  /**
122  * Get TDirectory pointer.
123  *
124  * The TFile pointer of an already opened file is recovered, else a new file is opened.
125  *
126  * \param id identifier
127  * \return pointer to TDirectory
128  */
129  inline TDirectory* getDirectory(const JRootObjectID& id)
130  {
131  TFile* in = getFile(id.getFilename().c_str(), "exist");
132 
133  if (in == NULL || !in->IsOpen()) {
134  return NULL;
135  }
136 
137  if (id.getDirectory() != "")
138  return in->GetDirectory(id.getDirectory());
139  else
140  return in;
141  }
142 
143 
144  /**
145  * Get first TObject with given identifier.
146  *
147  * \param id identifier
148  * \return pointer to TObject (or NULL)
149  */
150  inline TObject* getObject(const JRootObjectID& id)
151  {
152  TDirectory* dir = getDirectory(id);
153 
154  if (dir != NULL) {
155 
156  const TRegexp regexp(id.getObjectName());
157 
158  TIter iter(dir->GetListOfKeys());
159 
160  for (TKey* key; (key = (TKey*) iter.Next()) != NULL; ) {
161 
162  const TString tag(key->GetName());
163 
164  // option match
165 
166  if (tag.Index(regexp) != -1) {
167  return key->ReadObj();
168  }
169  }
170  }
171 
172  return NULL;
173  }
174 
175 
176  /**
177  * Get drawing option of TH1.
178  *
179  * \param object pointer to TObject
180  * \return true if TH1 looks like a line; else false
181  */
182  inline bool isTAttLine(const TObject* object)
183  {
184  const TH1* h1 = dynamic_cast<const TH1*>(object);
185 
186  if (h1 != NULL) {
187 
188  if (h1->GetSumw2N()) {
189 
190  for (Int_t i = 1; i <= h1->GetNbinsX(); ++i) {
191 
192  if (h1->GetBinError(i) != 0.0) {
193  return false;
194  }
195  }
196  }
197 
198  return true;
199 
200  } else {
201 
202  return false;
203  }
204  }
205 
206 
207  /**
208  * Get result of given textual formula.
209  *
210  * The formula may contain names of member methods of the object pointed to.\n
211  * These methods could have arguments and the return type should be compatible with <tt>Double_t</tt>.\n
212  * Example:
213  * <pre>
214  * getResult("1.0/GetEntries", TH1*);
215  * </pre>
216  *
217  * \param text text
218  * \param object pointer to object
219  * \return value
220  */
221  inline Double_t getResult(const TString& text, TObject* object = NULL)
222  {
223  TString buffer(text);
224 
225  if (object != NULL) {
226 
227  TClass* p = TClass::GetClass(object->ClassName());
228 
229  if (p != NULL) {
230 
231  for ( ; ; ) {
232 
233  TIterator* iter = p->GetListOfAllPublicMethods()->MakeIterator();
234  TMethod* method = NULL;
235 
236  for (TMethod* p; (p = (TMethod*) iter->Next()) != NULL; ) {
237  if (buffer.Index(p->GetName()) != -1) {
238  if (method == NULL || strlen(p->GetName()) > strlen(method->GetName())) {
239  method = p;
240  }
241  }
242  }
243 
244  if (method == NULL) {
245  break;
246  }
247 
248  for (Ssiz_t index; (index = buffer.Index(method->GetName())) != -1; ) {
249 
250  const TRegexp fp(" *([^)]*)"); // function call
251 
252  Ssiz_t len;
253  Ssiz_t pos = buffer.Index(fp, &len, index);
254 
255  Double_t value;
256 
257  if (pos == -1 || pos != index + (Ssiz_t) strlen(method->GetName())) {
258 
259  TMethodCall(p, method->GetName(), NULL).Execute(object, value);
260 
261  len = strlen(method->GetName());
262 
263  } else {
264 
265  TMethodCall(p, method->GetName(), NULL).Execute(object, TString(buffer(pos + 1, len - 2)), value);
266 
267  len += strlen(method->GetName());
268  }
269 
270  buffer.Replace(index, len, TString::Format("%15.5e", value));
271  }
272  }
273  }
274  }
275 
276  return TFormula("/tmp", buffer.Data()).Eval(0.0);
277  }
278 
279 
280  /**
281  * Get result of given textual formula.
282  *
283  * The formula may contain names of member methods of the object pointed to.\n
284  * These methods could have arguments and the return type should be compatible with <tt>Double_t</tt>.\n
285  * Example:
286  * <pre>
287  * getResult("1.0/GetEntries", TH1*);
288  * </pre>
289  *
290  * \param text text
291  * \param object pointer to object
292  * \return value
293  */
294  inline Double_t getResult(const std::string& text, TObject* object = NULL)
295  {
296  return getResult(TString(text.c_str()), object);
297  }
298 
299 
300  /**
301  * Get parameter number from text string.
302  *
303  * The number corresponds to the value <tt>[0-9]*</tt> in the expression <tt>"p[0-9]* = .."</tt>.
304  *
305  * \param text text
306  * \return parameter number
307  */
308  inline int getParameter(const std::string& text)
309  {
310  const char* regexp("p[0-9]* *=");
311 
312  TString buffer(text.c_str());
313 
314  buffer = buffer(TRegexp(regexp));
315  buffer = buffer(1, buffer.Length() - 2);
316 
317  if (!buffer.IsDigit()) {
318  THROW(JParseError, "Text is not a number " << text << ' ' << regexp);
319  }
320 
321  return buffer.Atoi();
322  }
323 
324 
325  /**
326  * Get parameter value from text string.
327  *
328  * The formula may contain names of member methods of the object pointed to.\n
329  * These methods could have arguments and the return type should be compatible with <tt>Double_t</tt>.\n
330  * Example:
331  * <pre>
332  * getValue("p[..] = 2*GetMaximum", TH1*);
333  * </pre>
334  *
335  * \param text text
336  * \param object pointer to object
337  * \return value
338  */
339  inline Double_t getValue(const std::string& text, TObject* object = NULL)
340  {
341  const char* regexp("=.*");
342 
343  TString buffer(text.c_str());
344 
345  buffer = buffer(TRegexp(regexp));
346  buffer = buffer(1, buffer.Length() - 1);
347 
348  return getResult(std::string(buffer), object);
349  }
350 
351 
352  /**
353  * Make histogram axis logarithmic (e.g.\ after filling with log10()).
354  *
355  * \param axis axis
356  */
357  inline void setLogarithmic(TAxis* axis)
358  {
359  if (axis != NULL) {
360 
361  const Int_t first = axis->GetFirst();
362  const Int_t last = axis->GetLast();
363 
364  const Double_t xmin = axis->GetBinLowEdge(first);
365  const Double_t xmax = axis->GetBinLowEdge(last) + axis->GetBinWidth(last);
366 
367  const Int_t N = axis->GetNbins();
368  Double_t buffer[N+1];
369 
370  buffer[0] = TMath::Power(10.0, axis->GetBinLowEdge(1));
371 
372  for (Int_t i = 1; i <= N; ++i) {
373  buffer[i] = TMath::Power(10.0, axis->GetBinLowEdge(i) + axis->GetBinWidth(i));
374  }
375 
376  axis->Set(N, buffer);
377 
378  if (axis->TestBit(TAxis::kAxisRange)) {
379  axis->SetRangeUser(TMath::Power(10.0, xmin), TMath::Power(10.0, xmax));
380  }
381  }
382  }
383 
384 
385  /**
386  * Make given parameter in formula logarithmic (e.g.\ after filling with log10()).
387  *
388  * \param formula formula
389  * \param parameter parameter
390  * \return formula
391  */
392  inline TString getLogarithmic(const TString& formula, const char parameter)
393  {
394  const TRegexp regexp[] = {
395  TString("^") + TString(parameter) + TString("[^a-zA-Z_]"), // parameter at start of line
396  TString("[^a-zA-Z_]") + TString(parameter) + TString("[^a-zA-Z_]"), // parameter in middle of line
397  TString("[^a-zA-Z_]") + TString(parameter) + TString("$") // parameter at end of line
398  };
399 
400  const TString replacement = TString("log10(") + TString(parameter) + TString(")");
401 
402  TString buffer(formula);
403 
404  if (buffer.Length() == 1 && buffer[0] == parameter) {
405 
406  buffer = replacement;
407 
408  } else {
409 
410  for (Ssiz_t pos = 0, i; pos < buffer.Length(); pos += replacement.Length()) {
411  if ((i = buffer.Index(regexp[0], pos)) != -1) { buffer.Replace((pos = i + 0), 1, replacement); }
412  else if ((i = buffer.Index(regexp[1], pos)) != -1) { buffer.Replace((pos = i + 1), 1, replacement); }
413  else if ((i = buffer.Index(regexp[2], pos)) != -1) { buffer.Replace((pos = i + 1), 1, replacement); }
414  else { break; }
415  }
416  }
417  return buffer;
418  }
419 
420 
421  /**
422  * Copy function parameters.
423  *
424  * \param from function
425  * \param to function
426  */
427  inline void copy(const TF1& from, TF1& to)
428  {
429  static_cast<TAttLine&> (to) = static_cast<const TAttLine&> (from);
430  static_cast<TAttFill&> (to) = static_cast<const TAttFill&> (from);
431  static_cast<TAttMarker&>(to) = static_cast<const TAttMarker&>(from);
432 
433  to.SetParameters(from.GetParameters());
434 
435  to.SetNpx(from.GetNpx());
436  }
437 
438 
439  /**
440  * Copy function parameters.
441  *
442  * \param from function
443  * \param to function
444  */
445  inline void copy(const TF2& from, TF2& to)
446  {
447  copy(static_cast<const TF1&>(from), static_cast<TF1&>(to));
448 
449  to.SetNpy(from.GetNpy());
450  }
451 
452 
453  /**
454  * Make parameter <tt>x</tt> of function logarithmic (e.g.\ after filling with log10()).
455  *
456  * \param f1 function
457  */
458  inline void setLogarithmicX(TF1* f1)
459  {
460  if (f1 != NULL) {
461 
462  TF1 fn(f1->GetName(), getLogarithmic(f1->GetExpFormula(), 'x'));
463 
464  copy(*f1, fn);
465 
466  fn.SetRange(f1->GetXmin(),
467  f1->GetXmax());
468 
469  *f1 = fn;
470 
471  f1->SetRange(TMath::Power(10.0,f1->GetXmin()),
472  TMath::Power(10.0,f1->GetXmax()));
473  }
474  }
475 
476 
477  /**
478  * Make parameter <tt>x</tt> of function logarithmic (e.g.\ after filling with log10()).
479  *
480  * \param f2 function
481  */
482  inline void setLogarithmicX(TF2* f2)
483  {
484  if (f2 != NULL) {
485 
486  TF2 fn(f2->GetName(), getLogarithmic(f2->GetExpFormula(), 'x'));
487 
488  copy(*f2, fn);
489 
490  fn.SetRange(f2->GetXmin(),
491  f2->GetYmin(),
492  f2->GetXmax(),
493  f2->GetYmax());
494 
495  *f2 = fn;
496 
497  f2->SetRange(TMath::Power(10.0,f2->GetXmin()),
498  f2->GetYmin(),
499  TMath::Power(10.0,f2->GetXmax()),
500  f2->GetYmax());
501  }
502  }
503 
504 
505  /**
506  * Make parameter <tt>y</tt> of function logarithmic (e.g.\ after filling with log10()).
507  *
508  * \param f2 function
509  */
510  inline void setLogarithmicY(TF2* f2)
511  {
512  if (f2 != NULL) {
513 
514  TF2 fn(f2->GetName(), getLogarithmic(f2->GetExpFormula(), 'y'));
515 
516  copy(*f2, fn);
517 
518  fn.SetRange(f2->GetXmin(),
519  f2->GetYmin(),
520  f2->GetXmax(),
521  f2->GetYmax());
522 
523  *f2 = fn;
524 
525  f2->SetRange(f2->GetXmin(),
526  TMath::Power(10.0,f2->GetYmin()),
527  f2->GetXmax(),
528  TMath::Power(10.0,f2->GetYmax()));
529  }
530  }
531 
532 
533  /**
534  * Make X-axis and associated functions of given histogram logarithmic (e.g.\ after filling with log10()).
535  *
536  * \param h1 histogram
537  */
538  inline void setLogarithmicX(TH1* h1)
539  {
540  if (h1 != NULL) {
541 
542  for (TIter iter(h1->GetListOfFunctions()); TF1* f1 = dynamic_cast<TF1*>(iter.Next()); ) {
543  setLogarithmicX(f1);
544  }
545 
546  setLogarithmic(h1->GetXaxis());
547  }
548  }
549 
550 
551  /**
552  * Make X-axis and associated functions of given histogram logarithmic (e.g.\ after filling with log10()).
553  *
554  * \param h2 histogram
555  */
556  inline void setLogarithmicX(TH2* h2)
557  {
558  using namespace std;
559 
560  if (h2 != NULL) {
561 
562  for (TIter iter(h2->GetListOfFunctions()); TF2* f2 = dynamic_cast<TF2*>(iter.Next()); ) {
563  setLogarithmicX(f2);
564  }
565 
566  setLogarithmic(h2->GetXaxis());
567  }
568  }
569 
570 
571  /**
572  * Make Y-axis and associated functions of given histogram logarithmic (e.g.\ after filling with log10()).
573  *
574  * \param h2 histogram
575  */
576  inline void setLogarithmicY(TH2* h2)
577  {
578  if (h2 != NULL) {
579 
580  for (TIter iter(h2->GetListOfFunctions()); TF2* f2 = dynamic_cast<TF2*>(iter.Next()); ) {
581  setLogarithmicY(f2);
582  }
583 
584  setLogarithmic(h2->GetYaxis());
585  }
586  }
587 
588 
589  /**
590  * Make parameter <tt>x</tt> of graph logarithmic (e.g.\ after filling with log10()).
591  *
592  * \param g1 graph
593  */
594  inline void setLogarithmicX(TGraph* g1)
595  {
596  if (g1 != NULL) {
597 
598  for (TIter iter(g1->GetListOfFunctions()); TF1* f1 = dynamic_cast<TF1*>(iter.Next()); ) {
599  setLogarithmicX(f1);
600  }
601 
602  for (Int_t i = 0; i != g1->GetN(); ++i) {
603  g1->GetX()[i] = pow(10.0, g1->GetX()[i]);
604  }
605  }
606  }
607 
608 
609  /**
610  * Make parameter <tt>x</tt> of graph logarithmic (e.g.\ after filling with log10()).
611  *
612  * \param g2 graph
613  */
614  inline void setLogarithmicX(TGraph2D* g2)
615  {
616  if (g2 != NULL) {
617 
618  for (TIter iter(g2->GetListOfFunctions()); TF2* f2 = dynamic_cast<TF2*>(iter.Next()); ) {
619  setLogarithmicX(f2);
620  }
621 
622  for (Int_t i = 0; i != g2->GetN(); ++i) {
623  g2->GetX()[i] = pow(10.0, g2->GetX()[i]);
624  }
625  }
626  }
627 
628 
629  /**
630  * Make parameter <tt>y</tt> of graph logarithmic (e.g.\ after filling with log10()).
631  *
632  * \param g2 graph
633  */
634  inline void setLogarithmicY(TGraph2D* g2)
635  {
636  if (g2 != NULL) {
637 
638  for (TIter iter(g2->GetListOfFunctions()); TF2* f2 = dynamic_cast<TF2*>(iter.Next()); ) {
639  setLogarithmicY(f2);
640  }
641 
642  for (Int_t i = 0; i != g2->GetN(); ++i) {
643  g2->GetY()[i] = pow(10.0, g2->GetY()[i]);
644  }
645  }
646  }
647 
648 
649  /**
650  * Make parameter <tt>x</tt> of multi graph logarithmic (e.g.\ after filling with log10()).
651  *
652  * \param m1 multi graph
653  */
654  inline void setLogarithmicX(TMultiGraph* m1)
655  {
656  if (m1 != NULL) {
657 
658  for (TIter i1(m1->GetListOfGraphs()); TGraph* g1 = dynamic_cast<TGraph*>(i1()); ) {
659  setLogarithmicX(g1);
660  }
661  }
662  }
663 
664 
665  /**
666  * Convert 1D histogram to PDF.
667  *
668  * Possible options are:
669  * - N normalise to histogram contents;
670  * - W divide by bin width;
671  * - E convert also bin errors.
672  *
673  * \param h1 histogram
674  * \param option option
675  * \param factor scaling factor
676  */
677  inline void convertToPDF(TH1& h1, const std::string& option = "NW", const double factor = 1.0)
678  {
679  using namespace std;
680 
681  const bool normalise = (option.find('N') != string::npos || option.find('n') != string::npos);
682  const bool bin_width = (option.find('W') != string::npos || option.find('w') != string::npos);
683  const bool use_error = (option.find('E') != string::npos || option.find('e') != string::npos);
684 
685  Double_t W = 1.0;
686 
687  if (normalise) {
688 
689  W = 0.0;
690 
691  for (Int_t i = 1; i <= h1.GetXaxis()->GetNbins(); ++i) {
692  W += h1.GetBinContent(i);
693  }
694  }
695 
696  if (W != 0.0) {
697 
698  for (Int_t i = 1; i <= h1.GetXaxis()->GetNbins(); ++i) {
699 
700  const Double_t w = W * (bin_width ? h1.GetXaxis()->GetBinWidth(i) : 1.0);
701 
702  h1.SetBinContent(i, h1.GetBinContent(i) * factor / w);
703 
704  if (use_error) {
705  h1.SetBinError(i, h1.GetBinError(i) * factor / w);
706  }
707  }
708  }
709  }
710 
711 
712  /**
713  * Convert 2D histogram to PDF.
714  *
715  * Possible options are:
716  * - N normalise to histogram contents;
717  * - X convert x-axis to PDF;
718  * - Y convert y-axis to PDF;
719  * - W divide by bin width;
720  * - E convert also bin errors.
721  *
722  * \param h2 histogram
723  * \param option option
724  * \param factor scaling factor
725  */
726  inline void convertToPDF(TH2& h2, const std::string& option = "NXYW", const double factor = 1.0)
727  {
728  using namespace std;
729 
730  const bool normalise = (option.find('N') != string::npos || option.find('n') != string::npos);
731  const bool X = (option.find('X') != string::npos || option.find('x') != string::npos);
732  const bool Y = (option.find('Y') != string::npos || option.find('y') != string::npos);
733  const bool bin_width = (option.find('W') != string::npos || option.find('w') != string::npos);
734  const bool use_error = (option.find('E') != string::npos || option.find('e') != string::npos);
735 
736  Double_t W = 1.0;
737 
738  if (X && Y) {
739 
740  if (normalise) {
741 
742  W = 0.0;
743 
744  for (Int_t ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
745  for (Int_t iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
746  W += h2.GetBinContent(ix,iy);
747  }
748  }
749  }
750 
751  if (W != 0.0) {
752 
753  for (Int_t ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
754  for (Int_t iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
755 
756  const Double_t w = W * (bin_width ? h2.GetXaxis()->GetBinWidth(ix) * h2.GetYaxis()->GetBinWidth(iy) : 1.0);
757 
758  h2.SetBinContent(ix, iy, h2.GetBinContent(ix,iy) * factor / w);
759 
760  if (use_error) {
761  h2.SetBinError(ix, iy, h2.GetBinError(ix,iy) * factor / w);
762  }
763  }
764  }
765  }
766 
767  } else if (X) {
768 
769  for (Int_t iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
770 
771  if (normalise) {
772 
773  W = 0.0;
774 
775  for (Int_t ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
776  W += h2.GetBinContent(ix,iy);
777  }
778  }
779 
780  if (W != 0.0) {
781 
782  for (Int_t ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
783 
784  const Double_t w = W * (bin_width ? h2.GetXaxis()->GetBinWidth(ix) : 1.0);
785 
786  h2.SetBinContent(ix, iy, h2.GetBinContent(ix,iy) * factor / w);
787 
788  if (use_error) {
789  h2.SetBinError(ix, iy, h2.GetBinError(ix,iy) * factor / w);
790  }
791  }
792  }
793  }
794 
795  } else if (Y) {
796 
797  for (Int_t ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
798 
799  if (normalise) {
800 
801  W = 0.0;
802 
803  for (Int_t iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
804  W += h2.GetBinContent(ix,iy);
805  }
806  }
807 
808  if (W != 0.0) {
809 
810  for (Int_t iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
811 
812  const Double_t w = W * (bin_width ? h2.GetYaxis()->GetBinWidth(iy) : 1.0);
813 
814  h2.SetBinContent(ix, iy, h2.GetBinContent(ix,iy) / w);
815 
816  if (use_error) {
817  h2.SetBinError(ix, iy, h2.GetBinError(ix,iy) / w);
818  }
819  }
820  }
821  }
822  }
823  }
824 
825 
826  /**
827  * Set limits of TGraph.
828  *
829  * \param g1 graph
830  */
831  inline void setLimits(TGraph& g1)
832  {
833  using namespace std;
834 
835  Double_t ymin = numeric_limits<Double_t>::max();
836  Double_t ymax = numeric_limits<Double_t>::lowest();
837 
838  for (Int_t i = 0; i != g1.GetN(); ++i) {
839 
840  const Double_t y = g1.GetY()[i];
841 
842  if (y > ymax) { ymax = y; }
843  if (y < ymin) { ymin = y; }
844  }
845 
846  g1.SetMinimum(ymin);
847  g1.SetMaximum(ymax);
848  }
849 
850 
851  /**
852  * Set limits of TGraph2D.
853  *
854  * \param g2 graph
855  */
856  inline void setLimits(TGraph2D& g2)
857  {
858  using namespace std;
859 
860  Double_t zmin = numeric_limits<Double_t>::max();
861  Double_t zmax = numeric_limits<Double_t>::lowest();
862 
863  for (Int_t i = 0; i != g2.GetN(); ++i) {
864 
865  const Double_t z = g2.GetZ()[i];
866 
867  if (z > zmax) { zmax = z; }
868  if (z < zmin) { zmin = z; }
869  }
870 
871  g2.SetMinimum(zmin);
872  g2.SetMaximum(zmax);
873  }
874 
875 
876  /**
877  * Set axis range.
878  *
879  * \param xmin lower limit (I/O)
880  * \param xmax upper limit (I/O)
881  * \param logx logarithmic
882  */
883  inline void setRange(double& xmin,
884  double& xmax,
885  const bool logx)
886  {
887  if (logx) {
888  xmin = log(xmin);
889  xmax = log(xmax);
890  }
891 
892  double dx = (xmax - xmin) * 0.1;
893 
894  if (xmin > dx || xmin < 0.0)
895  xmin -= dx;
896  else
897  xmin = 0.0;
898 
899  xmax += dx;
900 
901  if (logx) {
902  xmin = exp(xmin);
903  xmax = exp(xmax);
904  }
905  }
906 
907 
908  /**
909  * Set axis with PMT address labels.
910  *
911  * This methods sets the labels of the given axis to the sorted values of the PMT ring and position.\n
912  * It should normally be called before filling of the corresponding histogram.\n
913  * The filling should then be made with the textual representation of the PMT ring and position (i.e.\ JDETECTOR::JPMTPhysicalAddress::toString).
914  *
915  * Alternatively, the filling can be made with the index of the PMT in the address map of the corresponding module
916  * (i.e.\ JDETECTOR::JModuleAddressMap::getIndex).\n
917  * In that case, the method can be called before or after filling of the histogram.
918  *
919  * \param axis axis
920  * \param memo module address map
921  */
922  inline void setAxisLabels(TAxis *axis, const JModuleAddressMap& memo)
923  {
924  using namespace JPP;
925 
926  if (axis->GetNbins() == (int) memo.size()) {
927 
928  for (int i = 0; i != axis->GetNbins(); ++i) {
929 
930  const JPMTPhysicalAddress& address = memo[i];
931 
932  axis->SetBinLabel(i + 1, address.toString().c_str());
933  }
934 
935  } else {
936 
937  THROW(JValueOutOfRange, "Number of bins " << axis->GetNbins() << " != " << memo.size());
938  }
939  }
940 
941 
942  /**
943  * Set axis labels with PMT addresses.
944  *
945  * This methods sets the labels of the given axis to the sorted values of the PMT ring and position.\n
946  * It should be called after filling of the corresponding histogram.\n
947  * The filling should have been made with the PMT number (e.g.\ KM3NETDAQ::JDAQHit::getPMT).
948  *
949  * \param h1 histogram
950  * \param axis axis <tt>(x, X, y, Y, z, Z)</tt>
951  * \param memo module address map
952  */
953  inline void setAxisLabels(TH1& h1, const std::string axis, const JModuleAddressMap& memo)
954  {
955  using namespace JPP;
956 
957  TAxis* pax = NULL;
958 
959  if (axis == "x" || axis == "X")
960  pax = h1.GetXaxis();
961  else if (axis == "y" || axis == "Y")
962  pax = h1.GetYaxis();
963  else if (axis == "z" || axis == "Z")
964  pax = h1.GetZaxis();
965  else
966  THROW(JValueOutOfRange, "Invalid axis " << axis);
967 
968  if (pax->GetNbins() == (int) memo.size()) {
969 
970  for (int i = 0; i != pax->GetNbins(); ++i) {
971 
972  const JPMTPhysicalAddress& address = memo.getAddressTranslator(i);
973 
974  pax->SetBinLabel(i + 1, address.toString().c_str());
975  }
976 
977  h1.LabelsOption("a", axis.c_str()); // sort labels
978 
979  } else {
980 
981  THROW(JValueOutOfRange, "Number of bins " << pax->GetNbins() << " != " << memo.size());
982  }
983  }
984 }
985 
986 #endif
static const char *const sqrt()
Set contents to signed difference between squares.
data_type w[N+1][M+1]
Definition: JPolint.hh:741
Exceptions.
TObject * getObject(const JRootObjectID &id)
Get first TObject with given identifier.
double getValue(const JScale_t scale)
Get numerical value corresponding to scale.
Definition: JScale.hh:47
int getParameter(const std::string &text)
Get parameter number from text string.
static const char *const Divide()
ROOT TH1::Divide.
static const char *const Subtract()
ROOT TH1::Subtract.
static const char *const SAME_AS_INPUT()
Set name of output histogram to name of input histogram.
char text[TEXT_SIZE]
Definition: elog.cc:72
void setLogarithmicX(TF1 *f1)
Make parameter x of function logarithmic (e.g. after filling with log10()).
#define THROW(JException_t, A)
Marco for throwing exception with std::ostream compatible message.
Definition: JException.hh:670
Definition: JRoot.hh:19
static const char *const Replace()
Set contents to associated function.
void setLimits(TGraph &g1)
Set limits of TGraph.
static const char *const subtract()
Subtract contents with lookup bin in second histogram.
then for HISTOGRAM in h0 h1
Definition: JMatrixNZ.sh:69
Auxiliary class to handle file name, ROOT directory and object name.
then usage $script< string identifier >< detectorfile > input file(toashort file)+" "\nNote that the input files and toashort files should be one-to-one related." fi if (( $
then fatal Wrong number of arguments fi set_variable STRING $argv[1] set_variable DETECTORXY_TXT $WORKDIR $DETECTORXY_TXT tail read X Y CHI2 RMS printf optimum n $X $Y $CHI2 $RMS awk v Y
static const char *const stdev()
Set contents to standard deviation.
then set_variable FORMULA *[0] exp(-0.5 *(x-[1])*(x-[1])/([2]*[2]))*exp(-0.5 *(y-[1])*(y-[1])/([2]*[2]))" JF2 -o $WORKDIR/f2.root -F "$FORMULA" -@ "p0
then echo The file $DIR KM3NeT_00000001_00000000 root already please rename or remove it first
void setLogarithmic(TAxis *axis)
Make histogram axis logarithmic (e.g. after filling with log10()).
bool isTAttLine(const TObject *object)
Get drawing option of TH1.
static const char *const efficiency()
Divide contents and multiply errors with inefficiency.
static const char *const Add()
ROOT TH1::Add.
Lookup table for PMT addresses in optical module.
static const char *const multiply()
Multiply contents with lookup bin in second histogram.
TFile * getFile(const std::string &file_name, const std::string &option="exist")
Get TFile pointer corresponding to give file name.
static const char *const add()
Add contents with lookup bin in second histogram.
Auxiliary data structure for JOpera1D.cc and JOpera2D.cc applications.
Double_t getResult(const TString &text, TObject *object=NULL)
Get result of given textual formula.
void convertToPDF(TH1 &h1, const std::string &option="NW", const double factor=1.0)
Convert 1D histogram to PDF.
void setAxisLabels(TAxis *axis, const JModuleAddressMap &memo)
Set axis with PMT address labels.
T pow(const T &x, const double y)
Power .
Definition: JMath.hh:98
static const char *const Multiply()
ROOT TH1::Multiply.
static const char *const SAME_AS_OPERATION()
Set name of output histogram to name of operation.
void setRange(double &xmin, double &xmax, const bool logx)
Set axis range.
const JPMTAddressTranslator & getAddressTranslator(const int tdc) const
Get PMT address translator.
static const char *const divide()
Divide contents with lookup bin in second histogram.
void copy(const Head &from, JHead &to)
Copy header from from to to.
Definition: JHead.cc:139
Exception for parsing value.
Definition: JException.hh:180
TDirectory * getDirectory(const JRootObjectID &id)
Get TDirectory pointer.
std::string getFilename(const std::string &file_name)
Get file name part, i.e. part after last JEEP::PATHNAME_SEPARATOR if any.
Definition: JeepToolkit.hh:88
void setLogarithmicY(TF2 *f2)
Make parameter y of function logarithmic (e.g. after filling with log10()).
Exception for accessing a value in a collection that is outside of its range.
Definition: JException.hh:162
TString getLogarithmic(const TString &formula, const char parameter)
Make given parameter in formula logarithmic (e.g. after filling with log10()).
then fatal Wrong number of arguments fi set_variable DETECTOR $argv[1] set_variable INPUT_FILE $argv[2] eval JPrintDetector a $DETECTOR O IDENTIFIER eval JPrintDetector a $DETECTOR O SUMMARY source JAcoustics sh $DETECTOR_ID CHECK_EXIT_CODE typeset A TRIPODS get_tripods $WORKDIR tripod txt TRIPODS for EMITTER in
Definition: JCanberra.sh:38
static const char *const TIMESTAMP
Time stamp of earliest UTC time.
then usage $script[input file[working directory[option]]] nWhere option can be N
Definition: JMuonPostfit.sh:37
Double_t g1(const Double_t x)
Function.
Definition: JQuantiles.cc:25