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