Jpp  18.4.0
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 #include <memory>
9 
10 #pragma GCC diagnostic push
11 #pragma GCC diagnostic ignored "-Wall"
12 #include "TError.h"
13 #include "TFile.h"
14 #include "TClass.h"
15 #include "TObject.h"
16 #include "TKey.h"
17 #include "TH1.h"
18 #include "TH2.h"
19 #include "TH3.h"
20 #include "TGraph.h"
21 #include "TGraphErrors.h"
22 #include "TGraph2D.h"
23 #include "TGraph2DErrors.h"
24 #include "TMultiGraph.h"
25 #include "TLine.h"
26 #include "TEllipse.h"
27 #include "TString.h"
28 #include "TRegexp.h"
29 #include "TFormula.h"
30 #include "TF1.h"
31 #include "TF2.h"
32 #include "TIterator.h"
33 #include "TMethod.h"
34 #include "TMethodCall.h"
35 #include "TAxis.h"
36 #include "TMath.h"
37 #pragma GCC diagnostic pop
38 
39 #include "JLang/JException.hh"
40 #include "JGizmo/JRootObjectID.hh"
41 
43 
44 
45 /**
46  * \author mdejong
47  */
48 
49 namespace JGIZMO {}
50 namespace JPP { using namespace JGIZMO; }
51 
52 /**
53  * Auxiliary applications for use of ROOT and more.
54  */
55 namespace JGIZMO {
56 
57  using JLANG::JParseError;
61 
62 
63  /**
64  * Auxiliary data structure for JOpera1D.cc and JOpera2D.cc applications.
65  */
66  struct JOpera {
67  //
68  // Histogram name.
69  //
70  static const char* const SAME_AS_OPERATION() { return "%"; } //!< Set name of output histogram to name of operation
71  static const char* const SAME_AS_INPUT() { return "="; } //!< Set name of output histogram to name of input histogram
72 
73  //
74  // Histogram operations.
75  //
76  static const char* const Add() { return "Add"; } //!< ROOT TH1::Add
77  static const char* const add() { return "add"; } //!< Add contents with lookup bin in second histogram
78  static const char* const Subtract() { return "Subtract"; } //!< ROOT TH1::Subtract
79  static const char* const subtract() { return "subtract"; } //!< Subtract contents with lookup bin in second histogram
80  static const char* const Multiply() { return "Multiply"; } //!< ROOT TH1::Multiply
81  static const char* const multiply() { return "multiply"; } //!< Multiply contents with lookup bin in second histogram
82  static const char* const Divide() { return "Divide"; } //!< ROOT TH1::Divide
83  static const char* const divide() { return "divide"; } //!< Divide contents with lookup bin in second histogram
84  static const char* const efficiency() { return "efficiency"; } //!< Divide contents and multiply errors with inefficiency
85  static const char* const stdev() { return "stdev"; } //!< Set contents to standard deviation
86  static const char* const sqrt() { return "sqrt"; } //!< Set contents to signed difference between squares
87  static const char* const Replace() { return "Replace"; } //!< Set contents to associated function
88  };
89 
90 
91  /**
92  * Time stamp of earliest UTC time.
93  */
94  static const char* const TIMESTAMP = "#splitline{}{#splitline{%d:%m:%y}{ %H:%M}}%F1970-01-01 00:00:00";
95 
96 
97  /**
98  * Get TFile pointer corresponding to give file name.
99  *
100  * The TFile pointer of an already opened file is recovered, else a new file is opened.\n
101  * Note that the closure of the opened files should be done by the caller of this method.
102  *
103  * \param file_name file name
104  * \param option TFile::Open option
105  * \return pointer to TFile
106  */
107  inline TFile* getFile(const std::string& file_name, const std::string& option = "exist")
108  {
109  using namespace std;
110 
111  gErrorIgnoreLevel = kError;
112 
113  static map<string, TFile*> zmap;
114 
115  map<string, TFile*>::iterator i = zmap.find(file_name);
116 
117  if (i == zmap.end() || i->second == NULL || !i->second->IsOpen()) {
118 
119  TFile* file = TFile::Open(file_name.c_str(), option.c_str());
120 
121  zmap[file_name] = file;
122 
123  return file;
124 
125  } else {
126 
127  return i->second;
128  }
129  }
130 
131 
132  /**
133  * Get TDirectory pointer.
134  *
135  * The TFile pointer of an already opened file is recovered, else a new file is opened.
136  *
137  * \param id identifier
138  * \return pointer to TDirectory
139  */
140  inline TDirectory* getDirectory(const JRootObjectID& id)
141  {
142  TFile* in = getFile(id.getFilename().c_str(), "exist");
143 
144  if (in == NULL || !in->IsOpen()) {
145  return NULL;
146  }
147 
148  if (id.getDirectory() != "")
149  return in->GetDirectory(id.getDirectory());
150  else
151  return in;
152  }
153 
154 
155  /**
156  * Get first TObject with given identifier.
157  *
158  * \param id identifier
159  * \return pointer to TObject (or NULL)
160  */
161  inline TObject* getObject(const JRootObjectID& id)
162  {
163  TDirectory* dir = getDirectory(id);
164 
165  if (dir != NULL) {
166 
167  const TRegexp regexp(id.getObjectName());
168 
169  TIter iter(dir->GetListOfKeys());
170 
171  for (TKey* key; (key = (TKey*) iter.Next()) != NULL; ) {
172 
173  const TString tag(key->GetName());
174 
175  // option match
176 
177  if (tag.Index(regexp) != -1) {
178  return key->ReadObj();
179  }
180  }
181  }
182 
183  return NULL;
184  }
185 
186 
187  /**
188  * Get drawing option of TH1.
189  *
190  * \param object pointer to TObject
191  * \return true if TH1 looks like a line; else false
192  */
193  inline bool isTAttLine(const TObject* object)
194  {
195  {
196  const TH1* h1 = dynamic_cast<const TH1*>(object);
197 
198  if (h1 != NULL) {
199 
200  if (h1->GetSumw2N()) {
201  for (Int_t i = 1; i <= h1->GetNbinsX(); ++i) {
202  if (h1->GetBinError(i) != 0.0) {
203  return false;
204  }
205  }
206  }
207  }
208  }
209  {
210  const TGraphErrors* g1 = dynamic_cast<const TGraphErrors*>(object);
211 
212  if (g1 != NULL) {
213 
214  for (Int_t i = 0; i != g1->GetN(); ++i) {
215  if (g1->GetEY()[i] != 0.0) {
216  return false;
217  }
218  }
219 
220  return g1->GetN() > 1;
221  }
222  }
223  {
224  const TGraph* g1 = dynamic_cast<const TGraph*>(object);
225 
226  if (g1 != NULL) {
227  return g1->GetN() > 1;
228  }
229  }
230 
231  return true;
232  }
233 
234 
235  /**
236  * Get result of given textual formula.
237  *
238  * The formula may contain names of member methods of the object pointed to.\n
239  * These methods should return a value that is compatible with <tt>Double_t</tt> and could have arguments.\n
240  * For example:
241  * <pre>
242  * getResult("1.0 / GetEntries", TH1*);
243  * </pre>
244  *
245  * \param text text
246  * \param object pointer to object
247  * \return value
248  */
249  inline Double_t getResult(const TString& text, TObject* object = NULL)
250  {
251  TString buffer(text);
252 
253  if (object != NULL) {
254 
255  TClass* p = TClass::GetClass(object->ClassName());
256 
257  if (p != NULL) {
258 
259  for ( ; ; ) {
260 
261  TMethod* method = NULL;
262 
263  for (std::unique_ptr<TIterator> iter(p->GetListOfAllPublicMethods()->MakeIterator()); TMethod* p = (TMethod*) iter->Next(); ) {
264  if (buffer.Index(p->GetName()) != -1) {
265  if (method == NULL || strlen(p->GetName()) > strlen(method->GetName())) {
266  method = p;
267  }
268  }
269  }
270 
271  if (method == NULL) {
272  break;
273  }
274 
275  for (Ssiz_t index; (index = buffer.Index(method->GetName())) != -1; ) {
276 
277  const TRegexp fp(" *([^)]*)"); // function call
278 
279  Ssiz_t len;
280  Ssiz_t pos = buffer.Index(fp, &len, index);
281 
282  Double_t value;
283 
284  if (pos == -1 || pos != index + (Ssiz_t) strlen(method->GetName())) {
285 
286  TMethodCall(p, method->GetName(), NULL).Execute(object, value);
287 
288  len = strlen(method->GetName());
289 
290  } else {
291 
292  TMethodCall(p, method->GetName(), NULL).Execute(object, TString(buffer(pos + 1, len - 2)), value);
293 
294  len += strlen(method->GetName());
295  }
296 
297  buffer.Replace(index, len, TString::Format("%20.10e", value));
298  }
299  }
300  }
301  }
302 
303  return TFormula("/tmp", buffer.Data()).Eval(0.0);
304  }
305 
306 
307  /**
308  * Get result of given textual formula.
309  *
310  * The formula may contain names of member methods of the object pointed to.\n
311  * These methods should return a value that is compatible with <tt>Double_t</tt> and could have arguments.\n
312  * For example:
313  * <pre>
314  * getResult("1.0 / GetEntries", TH1*);
315  * </pre>
316  *
317  * \param text text
318  * \param object pointer to object
319  * \return value
320  */
321  inline Double_t getResult(const std::string& text, TObject* object = NULL)
322  {
323  return getResult(TString(text.c_str()), object);
324  }
325 
326 
327  /**
328  * Get parameter number from text string.
329  *
330  * The number corresponds to the value <tt>[0-9]*</tt> in the expression <tt>"p[0-9]* = .."</tt>.
331  *
332  * \param text text
333  * \return parameter number
334  */
335  inline int getParameter(const std::string& text)
336  {
337  const char* regexp("p[0-9]* *=");
338 
339  TString buffer(text.c_str());
340 
341  buffer = buffer(TRegexp(regexp));
342  buffer = buffer(1, buffer.Length() - 2);
343 
344  if (!buffer.IsDigit()) {
345  THROW(JParseError, "Text is not a number " << text << ' ' << regexp);
346  }
347 
348  return buffer.Atoi();
349  }
350 
351 
352  /**
353  * Get parameter value from text string.
354  *
355  * The formula may contain names of member methods of the object pointed to.\n
356  * These methods should return a value that is compatible with <tt>Double_t</tt> and could have arguments.\n
357  * For example:
358  * <pre>
359  * getValue("p[..] = 2 * GetMaximum", TH1*);
360  * </pre>
361  *
362  * \param text text
363  * \param object pointer to object
364  * \return value
365  */
366  inline Double_t getValue(const std::string& text, TObject* object = NULL)
367  {
368  const char* regexp("=.*");
369 
370  TString buffer(text.c_str());
371 
372  buffer = buffer(TRegexp(regexp));
373  buffer = buffer(1, buffer.Length() - 1);
374 
375  return getResult(std::string(buffer), object);
376  }
377 
378 
379  /**
380  * Get parameter value from text string.
381  *
382  * The formula may contain names of member methods of the object pointed to.\n
383  * These methods should return a value that is compatible with <tt>Double_t</tt> and could have arguments.\n
384  * For example:
385  * <pre>
386  * getValue("p[..] = 1.0 2.0 3.0", 1);
387  * </pre>
388  * will return <tt>2.0</tt>.
389  *
390  * \param text text
391  * \param index index
392  * \return value
393  */
394  inline Double_t getValue(const std::string& text, const int index)
395  {
396  using namespace std;
397 
398  const char* regexp("=.*");
399 
400  TString buffer(text.c_str());
401 
402  buffer = buffer(TRegexp(regexp));
403  buffer = buffer(1, buffer.Length() - 1);
404 
405 
406  istringstream is((std::string(buffer)));
407 
408  Double_t value;
409 
410  for (int i = index; is >> value && i > 0; --i) {}
411 
412  if (is)
413  return value;
414  else
415  THROW(JParseError, "Text does not contain a number at given position " << buffer << ' ' << index);
416  }
417 
418 
419  /**
420  * Make histogram axis logarithmic (e.g.\ after using <tt>log10()</tt>).
421  *
422  * \param axis axis
423  */
424  inline void setLogarithmic(TAxis* axis)
425  {
426  if (axis != NULL) {
427 
428  const Int_t first = axis->GetFirst();
429  const Int_t last = axis->GetLast();
430 
431  const Double_t xmin = axis->GetBinLowEdge(first);
432  const Double_t xmax = axis->GetBinLowEdge(last) + axis->GetBinWidth(last);
433 
434  const Int_t N = axis->GetNbins();
435  Double_t buffer[N+1];
436 
437  buffer[0] = TMath::Power(10.0, axis->GetBinLowEdge(1));
438 
439  for (Int_t i = 1; i <= N; ++i) {
440  buffer[i] = TMath::Power(10.0, axis->GetBinLowEdge(i) + axis->GetBinWidth(i));
441  }
442 
443  axis->Set(N, buffer);
444 
445  if (axis->TestBit(TAxis::kAxisRange)) {
446  axis->SetRangeUser(TMath::Power(10.0, xmin), TMath::Power(10.0, xmax));
447  }
448  }
449  }
450 
451 
452  /**
453  * Make given parameter in formula logarithmic (e.g.\ after using <tt>log10()</tt>).
454  *
455  * \param formula formula
456  * \param parameter parameter
457  * \return formula
458  */
459  inline TString getLogarithmic(const TString& formula, const char parameter)
460  {
461  const TRegexp regexp[] = {
462  TString("^") + TString(parameter) + TString("[^a-zA-Z_]"), // parameter at start of line
463  TString("[^a-zA-Z_]") + TString(parameter) + TString("[^a-zA-Z_]"), // parameter in middle of line
464  TString("[^a-zA-Z_]") + TString(parameter) + TString("$") // parameter at end of line
465  };
466 
467  const TString replacement = TString("log10(") + TString(parameter) + TString(")");
468 
469  TString buffer(formula);
470 
471  if (buffer.Length() == 1 && buffer[0] == parameter) {
472 
473  buffer = replacement;
474 
475  } else {
476 
477  for (Ssiz_t pos = 0, i; pos < buffer.Length(); pos += replacement.Length()) {
478  if ((i = buffer.Index(regexp[0], pos)) != -1) { buffer.Replace((pos = i + 0), 1, replacement); }
479  else if ((i = buffer.Index(regexp[1], pos)) != -1) { buffer.Replace((pos = i + 1), 1, replacement); }
480  else if ((i = buffer.Index(regexp[2], pos)) != -1) { buffer.Replace((pos = i + 1), 1, replacement); }
481  else { break; }
482  }
483  }
484 
485  return buffer;
486  }
487 
488 
489  /**
490  * Copy function parameters.
491  *
492  * \param from function
493  * \param to function
494  */
495  inline void copy(const TF1& from, TF1& to)
496  {
497  static_cast<TAttLine&> (to) = static_cast<const TAttLine&> (from);
498  static_cast<TAttFill&> (to) = static_cast<const TAttFill&> (from);
499  static_cast<TAttMarker&>(to) = static_cast<const TAttMarker&>(from);
500 
501  to.SetParameters(from.GetParameters());
502 
503  to.SetNpx(from.GetNpx());
504  }
505 
506 
507  /**
508  * Copy function parameters.
509  *
510  * \param from function
511  * \param to function
512  */
513  inline void copy(const TF2& from, TF2& to)
514  {
515  copy(static_cast<const TF1&>(from), static_cast<TF1&>(to));
516 
517  to.SetNpy(from.GetNpy());
518  }
519 
520 
521  /**
522  * Make x-axis of objects in list logarithmic (e.g.\ after using <tt>log10()</tt>).
523  *
524  * \param list list
525  */
526  template<class T>
527  inline void setLogarithmicX(TList* list);
528 
529 
530  /**
531  * Make y-axis of objects in list logarithmic (e.g.\ after using <tt>log10()</tt>).
532  *
533  * \param list list
534  */
535  template<class T>
536  inline void setLogarithmicY(TList* list);
537 
538 
539  /**
540  * Make x-axis of given function logarithmic (e.g.\ after using <tt>log10()</tt>).
541  *
542  * \param f1 function
543  */
544  inline void setLogarithmicX(TF1* f1)
545  {
546  if (f1 != NULL) {
547 
548  TF1 fn(f1->GetName(), getLogarithmic(f1->GetExpFormula(), 'x'));
549 
550  copy(*f1, fn);
551 
552  fn.SetRange(f1->GetXmin(),
553  f1->GetXmax());
554 
555  *f1 = fn;
556 
557  f1->SetRange(TMath::Power(10.0,f1->GetXmin()),
558  TMath::Power(10.0,f1->GetXmax()));
559  }
560  }
561 
562 
563  /**
564  * Make y-axis of given function logarithmic (e.g.\ after using <tt>log10()</tt>).
565  *
566  * \param f1 function
567  */
568  inline void setLogarithmicY(TF1* f1)
569  {
570  if (f1 != NULL) {
571 
572  TString buffer = f1->GetExpFormula();
573 
574  buffer = "pow(10.0, " + buffer + ")";
575 
576  TF1 fn(f1->GetName(), buffer);
577 
578  copy(*f1, fn);
579 
580  *f1 = fn;
581  }
582  }
583 
584 
585  /**
586  * Make x-axis of given function logarithmic (e.g.\ after using <tt>log10()</tt>).
587  *
588  * \param f2 function
589  */
590  inline void setLogarithmicX(TF2* f2)
591  {
592  if (f2 != NULL) {
593 
594  TF2 fn(f2->GetName(), getLogarithmic(f2->GetExpFormula(), 'x'));
595 
596  copy(*f2, fn);
597 
598  fn.SetRange(f2->GetXmin(),
599  f2->GetYmin(),
600  f2->GetXmax(),
601  f2->GetYmax());
602 
603  *f2 = fn;
604 
605  f2->SetRange(TMath::Power(10.0,f2->GetXmin()),
606  f2->GetYmin(),
607  TMath::Power(10.0,f2->GetXmax()),
608  f2->GetYmax());
609  }
610  }
611 
612 
613  /**
614  * Make y-axis of given function logarithmic (e.g.\ after using <tt>log10()</tt>).
615  *
616  * \param f2 function
617  */
618  inline void setLogarithmicY(TF2* f2)
619  {
620  if (f2 != NULL) {
621 
622  TF2 fn(f2->GetName(), getLogarithmic(f2->GetExpFormula(), 'y'));
623 
624  copy(*f2, fn);
625 
626  fn.SetRange(f2->GetXmin(),
627  f2->GetYmin(),
628  f2->GetXmax(),
629  f2->GetYmax());
630 
631  *f2 = fn;
632 
633  f2->SetRange(f2->GetXmin(),
634  TMath::Power(10.0,f2->GetYmin()),
635  f2->GetXmax(),
636  TMath::Power(10.0,f2->GetYmax()));
637  }
638  }
639 
640 
641  /**
642  * Make x-axis and associated functions of given histogram logarithmic (e.g.\ after using <tt>log10()</tt>).
643  *
644  * \param h1 histogram
645  */
646  inline void setLogarithmicX(TH1* h1)
647  {
648  if (h1 != NULL) {
649 
650  setLogarithmicX<TF1>(h1->GetListOfFunctions());
651 
652  setLogarithmic(h1->GetXaxis());
653  }
654  }
655 
656 
657  /**
658  * Make x-axis and associated functions of given histogram logarithmic (e.g.\ after using <tt>log10()</tt>).
659  *
660  * \param h2 histogram
661  */
662  inline void setLogarithmicX(TH2* h2)
663  {
664  using namespace std;
665 
666  if (h2 != NULL) {
667 
668  setLogarithmicX<TF2>(h2->GetListOfFunctions());
669 
670  setLogarithmic(h2->GetXaxis());
671  }
672  }
673 
674 
675  /**
676  * Make y-axis and associated functions of given histogram logarithmic (e.g.\ after using <tt>log10()</tt>).
677  *
678  * \param h2 histogram
679  */
680  inline void setLogarithmicY(TH2* h2)
681  {
682  if (h2 != NULL) {
683 
684  setLogarithmicY<TF2>(h2->GetListOfFunctions());
685 
686  setLogarithmic(h2->GetYaxis());
687  }
688  }
689 
690 
691  /**
692  * Make x-axis and associated functions of given graph logarithmic (e.g.\ after using <tt>log10()</tt>).
693  *
694  * \param g1 graph
695  */
696  inline void setLogarithmicX(TGraph* g1)
697  {
698  if (g1 != NULL) {
699 
700  setLogarithmicX<TF1>(g1->GetListOfFunctions());
701 
702  for (Int_t i = 0; i != g1->GetN(); ++i) {
703  g1->GetX()[i] = pow(10.0, g1->GetX()[i]);
704  }
705  }
706  }
707 
708 
709  /**
710  * Make y-axis and associated functions of given graph logarithmic (e.g.\ after using <tt>log10()</tt>).
711  *
712  * \param g1 graph
713  */
714  inline void setLogarithmicY(TGraph* g1)
715  {
716  if (g1 != NULL) {
717 
718  setLogarithmicY<TF1>(g1->GetListOfFunctions());
719 
720  for (Int_t i = 0; i != g1->GetN(); ++i) {
721  g1->GetY()[i] = pow(10.0, g1->GetY()[i]);
722  }
723  }
724  }
725 
726 
727  /**
728  * Make x-axis and associated functions of given graph logarithmic (e.g.\ after using <tt>log10()</tt>).
729  *
730  * \param g1 graph
731  */
732  inline void setLogarithmicX(TGraphErrors* g1)
733  {
734  if (g1 != NULL) {
735 
736  setLogarithmicX<TF1>(g1->GetListOfFunctions());
737 
738  for (Int_t i = 0; i != g1->GetN(); ++i) {
739  g1->GetEX()[i] = pow(10.0, g1->GetX()[i] + g1->GetEX()[i]) - pow(10.0, g1->GetX()[i]);
740  g1->GetX() [i] = pow(10.0, g1->GetX()[i]);
741  }
742  }
743  }
744 
745 
746  /**
747  * Make y-axis and associated functions of given graph logarithmic (e.g.\ after using <tt>log10()</tt>).
748  *
749  * \param g1 graph
750  */
751  inline void setLogarithmicY(TGraphErrors* g1)
752  {
753  if (g1 != NULL) {
754 
755  setLogarithmicY<TF1>(g1->GetListOfFunctions());
756 
757  for (Int_t i = 0; i != g1->GetN(); ++i) {
758  g1->GetEY()[i] = pow(10.0, g1->GetY()[i] + g1->GetEY()[i]) - pow(10.0, g1->GetY()[i]);
759  g1->GetY() [i] = pow(10.0, g1->GetY()[i]);
760  }
761  }
762  }
763 
764 
765  /**
766  * Make x-axis and associated functions of given graph logarithmic (e.g.\ after using <tt>log10()</tt>).
767  *
768  * \param g2 graph
769  */
770  inline void setLogarithmicX(TGraph2D* g2)
771  {
772  if (g2 != NULL) {
773 
774  setLogarithmicX<TF2>(g2->GetListOfFunctions());
775 
776  for (Int_t i = 0; i != g2->GetN(); ++i) {
777  g2->GetX()[i] = pow(10.0, g2->GetX()[i]);
778  }
779  }
780  }
781 
782 
783  /**
784  * Make y-axis and associated functions of given graph logarithmic (e.g.\ after using <tt>log10()</tt>).
785  *
786  * \param g2 graph
787  */
788  inline void setLogarithmicY(TGraph2D* g2)
789  {
790  if (g2 != NULL) {
791 
792  setLogarithmicY<TF2>(g2->GetListOfFunctions());
793 
794  for (Int_t i = 0; i != g2->GetN(); ++i) {
795  g2->GetY()[i] = pow(10.0, g2->GetY()[i]);
796  }
797  }
798  }
799 
800 
801  /**
802  * Make x-axis and associated functions of given graph logarithmic (e.g.\ after using <tt>log10()</tt>).
803  *
804  * \param g2 graph
805  */
806  inline void setLogarithmicX(TGraph2DErrors* g2)
807  {
808  if (g2 != NULL) {
809 
810  setLogarithmicX<TF2>(g2->GetListOfFunctions());
811 
812  for (Int_t i = 0; i != g2->GetN(); ++i) {
813  g2->GetEX()[i] = pow(10.0, g2->GetX()[i] + g2->GetEX()[i]) - pow(10.0, g2->GetX()[i]);
814  g2->GetX() [i] = pow(10.0, g2->GetX()[i]);
815  }
816  }
817  }
818 
819 
820  /**
821  * Make y-axis and associated functions of given graph logarithmic (e.g.\ after using <tt>log10()</tt>).
822  *
823  * \param g2 graph
824  */
825  inline void setLogarithmicY(TGraph2DErrors* g2)
826  {
827  if (g2 != NULL) {
828 
829  setLogarithmicY<TF2>(g2->GetListOfFunctions());
830 
831  for (Int_t i = 0; i != g2->GetN(); ++i) {
832  g2->GetEY()[i] = pow(10.0, g2->GetY()[i] + g2->GetEY()[i]) - pow(10.0, g2->GetY()[i]);
833  g2->GetY() [i] = pow(10.0, g2->GetY()[i]);
834  }
835  }
836  }
837 
838 
839  /**
840  * Make x-axis of given multi-graph logarithmic (e.g.\ after using <tt>log10()</tt>).
841  *
842  * \param gn multi graph
843  */
844  inline void setLogarithmicX(TMultiGraph* gn)
845  {
846  if (gn != NULL) {
847  setLogarithmicX<TGraph>(gn->GetListOfGraphs());
848  }
849  }
850 
851 
852  /**
853  * Make y-axis of given multi-graph logarithmic (e.g.\ after using <tt>log10()</tt>).
854  *
855  * \param gn multi graph
856  */
857  inline void setLogarithmicY(TMultiGraph* gn)
858  {
859  if (gn != NULL) {
860  setLogarithmicY<TGraph>(gn->GetListOfGraphs());
861  }
862  }
863 
864 
865  /**
866  * Make x-axis of given line logarithmic (e.g.\ after using <tt>log10()</tt>).
867  *
868  * \param line line
869  */
870  inline void setLogarithmicX(TLine* line)
871  {
872  if (line != NULL) {
873  line->SetX1(pow(10.0, line->GetX1()));
874  line->SetX2(pow(10.0, line->GetX2()));
875  }
876  }
877 
878 
879  /**
880  * Make y-axis of given line logarithmic (e.g.\ after using <tt>log10()</tt>).
881  *
882  * \param line line
883  */
884  inline void setLogarithmicY(TLine* line)
885  {
886  if (line != NULL) {
887  line->SetY1(pow(10.0, line->GetY1()));
888  line->SetY2(pow(10.0, line->GetY2()));
889  }
890  }
891 
892 
893  /**
894  * Make x-axis of given ellipse logarithmic (e.g.\ after using <tt>log10()</tt>).
895  *
896  * \param ellipse ellipse
897  */
898  inline void setLogarithmicX(TEllipse* ellipse)
899  {
900  THROW(JFunctionalException, "Operation setLogarithmicX on TEllipse not allowed.");
901  }
902 
903 
904  /**
905  * Make y-axis of given ellipse logarithmic (e.g.\ after using <tt>log10()</tt>).
906  *
907  * \param ellipse ellipse
908  */
909  inline void setLogarithmicY(TEllipse* ellipse)
910  {
911  THROW(JFunctionalException, "Operation setLogarithmicY on TEllipse not allowed.");
912  }
913 
914 
915  /**
916  * Make x-axis of objects in list logarithmic (e.g.\ after using <tt>log10()</tt>).
917  *
918  * \param list list
919  */
920  template<class T>
921  inline void setLogarithmicX(TList* list)
922  {
923  for (TIter i(list); T* p = dynamic_cast<T*>(i.Next()); ) {
924  setLogarithmicX(p);
925  }
926  }
927 
928 
929  /**
930  * Make y-axis of objects in list logarithmic (e.g.\ after using <tt>log10()</tt>).
931  *
932  * \param list list
933  */
934  template<class T>
935  inline void setLogarithmicY(TList* list)
936  {
937  for (TIter i(list); T* p = dynamic_cast<T*>(i.Next()); ) {
938  setLogarithmicY(p);
939  }
940  }
941 
942 
943  /**
944  * Convert 1D histogram to PDF.
945  *
946  * Possible options are:
947  * - N normalise to histogram contents;
948  * - W divide by bin width;
949  * - E convert also bin errors.
950  *
951  * \param h1 histogram
952  * \param option option
953  * \param factor scaling factor
954  */
955  inline void convertToPDF(TH1& h1, const std::string& option = "NW", const double factor = 1.0)
956  {
957  using namespace std;
958 
959  const bool normalise = (option.find('N') != string::npos || option.find('n') != string::npos);
960  const bool bin_width = (option.find('W') != string::npos || option.find('w') != string::npos);
961  const bool use_error = (option.find('E') != string::npos || option.find('e') != string::npos);
962 
963  Double_t W = 1.0;
964 
965  if (normalise) {
966 
967  W = 0.0;
968 
969  for (Int_t i = 1; i <= h1.GetXaxis()->GetNbins(); ++i) {
970  W += h1.GetBinContent(i);
971  }
972  }
973 
974  if (W != 0.0) {
975 
976  for (Int_t i = 1; i <= h1.GetXaxis()->GetNbins(); ++i) {
977 
978  const Double_t w = W * (bin_width ? h1.GetXaxis()->GetBinWidth(i) : 1.0);
979 
980  h1.SetBinContent(i, h1.GetBinContent(i) * factor / w);
981 
982  if (use_error) {
983  h1.SetBinError(i, h1.GetBinError(i) * factor / w);
984  }
985  }
986  }
987  }
988 
989 
990  /**
991  * Convert 2D histogram to PDF.
992  *
993  * Possible options are:
994  * - N normalise to histogram contents;
995  * - X convert x-axis to PDF;
996  * - Y convert y-axis to PDF;
997  * - W divide by bin width;
998  * - E convert also bin errors.
999  *
1000  * \param h2 histogram
1001  * \param option option
1002  * \param factor scaling factor
1003  */
1004  inline void convertToPDF(TH2& h2, const std::string& option = "NXYW", const double factor = 1.0)
1005  {
1006  using namespace std;
1007 
1008  const bool normalise = (option.find('N') != string::npos || option.find('n') != string::npos);
1009  const bool X = (option.find('X') != string::npos || option.find('x') != string::npos);
1010  const bool Y = (option.find('Y') != string::npos || option.find('y') != string::npos);
1011  const bool bin_width = (option.find('W') != string::npos || option.find('w') != string::npos);
1012  const bool use_error = (option.find('E') != string::npos || option.find('e') != string::npos);
1013 
1014  Double_t W = 1.0;
1015 
1016  if (X && Y) {
1017 
1018  if (normalise) {
1019 
1020  W = 0.0;
1021 
1022  for (Int_t ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
1023  for (Int_t iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
1024  W += h2.GetBinContent(ix,iy);
1025  }
1026  }
1027  }
1028 
1029  if (W != 0.0) {
1030 
1031  for (Int_t ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
1032  for (Int_t iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
1033 
1034  const Double_t w = W * (bin_width ? h2.GetXaxis()->GetBinWidth(ix) * h2.GetYaxis()->GetBinWidth(iy) : 1.0);
1035 
1036  h2.SetBinContent(ix, iy, h2.GetBinContent(ix,iy) * factor / w);
1037 
1038  if (use_error) {
1039  h2.SetBinError(ix, iy, h2.GetBinError(ix,iy) * factor / w);
1040  }
1041  }
1042  }
1043  }
1044 
1045  } else if (X) {
1046 
1047  for (Int_t iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
1048 
1049  if (normalise) {
1050 
1051  W = 0.0;
1052 
1053  for (Int_t ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
1054  W += h2.GetBinContent(ix,iy);
1055  }
1056  }
1057 
1058  if (W != 0.0) {
1059 
1060  for (Int_t ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
1061 
1062  const Double_t w = W * (bin_width ? h2.GetXaxis()->GetBinWidth(ix) : 1.0);
1063 
1064  h2.SetBinContent(ix, iy, h2.GetBinContent(ix,iy) * factor / w);
1065 
1066  if (use_error) {
1067  h2.SetBinError(ix, iy, h2.GetBinError(ix,iy) * factor / w);
1068  }
1069  }
1070  }
1071  }
1072 
1073  } else if (Y) {
1074 
1075  for (Int_t ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
1076 
1077  if (normalise) {
1078 
1079  W = 0.0;
1080 
1081  for (Int_t iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
1082  W += h2.GetBinContent(ix,iy);
1083  }
1084  }
1085 
1086  if (W != 0.0) {
1087 
1088  for (Int_t iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
1089 
1090  const Double_t w = W * (bin_width ? h2.GetYaxis()->GetBinWidth(iy) : 1.0);
1091 
1092  h2.SetBinContent(ix, iy, h2.GetBinContent(ix,iy) / w);
1093 
1094  if (use_error) {
1095  h2.SetBinError(ix, iy, h2.GetBinError(ix,iy) / w);
1096  }
1097  }
1098  }
1099  }
1100  }
1101  }
1102 
1103 
1104  /**
1105  * Convert 3D histogram to PDF.
1106  *
1107  * Possible options are:
1108  * - N normalise to histogram contents;
1109  * - X convert x-axis to PDF;
1110  * - Y convert y-axis to PDF;
1111  * - Z convert z-axis to PDF;
1112  * - W divide by bin width;
1113  * - E convert also bin errors.
1114  *
1115  * \param h3 histogram
1116  * \param option option
1117  * \param factor scaling factor
1118  */
1119  inline void convertToPDF(TH3& h3, const std::string& option = "NXYW", const double factor = 1.0)
1120  {
1121  using namespace std;
1122 
1123  const bool normalise = (option.find('N') != string::npos || option.find('n') != string::npos);
1124  const bool X = (option.find('X') != string::npos || option.find('x') != string::npos);
1125  const bool Y = (option.find('Y') != string::npos || option.find('y') != string::npos);
1126  const bool Z = (option.find('Z') != string::npos || option.find('z') != string::npos);
1127  const bool bin_width = (option.find('W') != string::npos || option.find('w') != string::npos);
1128  const bool use_error = (option.find('E') != string::npos || option.find('e') != string::npos);
1129 
1130  Double_t W = 1.0;
1131 
1132  if (X && Y && Z) {
1133 
1134  if (normalise) {
1135 
1136  W = 0.0;
1137 
1138  for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1139  for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1140  for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1141  W += h3.GetBinContent(ix,iy,iz);
1142  }
1143  }
1144  }
1145  }
1146 
1147  if (W != 0.0) {
1148 
1149  for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1150  for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1151  for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1152 
1153  const Double_t w = W * (bin_width ? h3.GetXaxis()->GetBinWidth(ix) * h3.GetYaxis()->GetBinWidth(iy) * h3.GetZaxis()->GetBinWidth(iz) : 1.0);
1154 
1155  h3.SetBinContent(ix, iy, iz, h3.GetBinContent(ix,iy,iz) * factor / w);
1156 
1157  if (use_error) {
1158  h3.SetBinError(ix, iy, iz, h3.GetBinError(ix,iy,iz) * factor / w);
1159  }
1160  }
1161  }
1162  }
1163  }
1164 
1165  } else if (X && Z) {
1166 
1167  for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1168 
1169  if (normalise) {
1170 
1171  W = 0.0;
1172 
1173  for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1174  for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1175  W += h3.GetBinContent(ix,iy,iz);
1176  }
1177  }
1178  }
1179 
1180  if (W != 0.0) {
1181 
1182  for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1183  for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1184 
1185  const Double_t w = W * (bin_width ? h3.GetXaxis()->GetBinWidth(ix) * h3.GetZaxis()->GetBinWidth(iz) : 1.0);
1186 
1187  h3.SetBinContent(ix, iy, iz, h3.GetBinContent(ix,iy,iz) * factor / w);
1188 
1189  if (use_error) {
1190  h3.SetBinError(ix, iy, iz, h3.GetBinError(ix,iy,iz) * factor / w);
1191  }
1192  }
1193  }
1194  }
1195  }
1196 
1197  } else if (Y && Z) {
1198 
1199  for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1200 
1201  if (normalise) {
1202 
1203  W = 0.0;
1204 
1205  for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1206  for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1207  W += h3.GetBinContent(ix,iy,iz);
1208  }
1209  }
1210  }
1211 
1212  if (W != 0.0) {
1213 
1214  for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1215  for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1216 
1217  const Double_t w = W * (bin_width ? h3.GetYaxis()->GetBinWidth(iy) * h3.GetZaxis()->GetBinWidth(iz) : 1.0);
1218 
1219  h3.SetBinContent(ix, iy, iz, h3.GetBinContent(ix,iy,iz) * factor / w);
1220 
1221  if (use_error) {
1222  h3.SetBinError(ix, iy, iz, h3.GetBinError(ix,iy,iz) * factor / w);
1223  }
1224  }
1225  }
1226  }
1227  }
1228 
1229  } else if (X && Y) {
1230 
1231  for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1232 
1233  if (normalise) {
1234 
1235  W = 0.0;
1236 
1237  for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1238  for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1239  W += h3.GetBinContent(ix,iy,iz);
1240  }
1241  }
1242  }
1243 
1244  if (W != 0.0) {
1245 
1246  for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1247  for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1248 
1249  const Double_t w = W * (bin_width ? h3.GetXaxis()->GetBinWidth(ix) * h3.GetYaxis()->GetBinWidth(iy) : 1.0);
1250 
1251  h3.SetBinContent(ix, iy, iz, h3.GetBinContent(ix,iy,iz) * factor / w);
1252 
1253  if (use_error) {
1254  h3.SetBinError(ix, iy, iz, h3.GetBinError(ix,iy,iz) * factor / w);
1255  }
1256  }
1257  }
1258  }
1259  }
1260 
1261  } else if (X) {
1262 
1263  for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1264  for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1265 
1266  if (normalise) {
1267 
1268  W = 0.0;
1269 
1270  for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1271  W += h3.GetBinContent(ix,iy,iz);
1272  }
1273  }
1274 
1275  if (W != 0.0) {
1276 
1277  for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1278 
1279  const Double_t w = W * (bin_width ? h3.GetXaxis()->GetBinWidth(ix) : 1.0);
1280 
1281  h3.SetBinContent(ix, iy, iz, h3.GetBinContent(ix,iy,iz) * factor / w);
1282 
1283  if (use_error) {
1284  h3.SetBinError(ix, iy, iz, h3.GetBinError(ix,iy,iz) * factor / w);
1285  }
1286  }
1287  }
1288  }
1289  }
1290 
1291  } else if (Y) {
1292 
1293  for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1294  for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1295 
1296  if (normalise) {
1297 
1298  W = 0.0;
1299 
1300  for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1301  W += h3.GetBinContent(ix,iy,iz);
1302  }
1303  }
1304 
1305  if (W != 0.0) {
1306 
1307  for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1308 
1309  const Double_t w = W * (bin_width ? h3.GetYaxis()->GetBinWidth(iy) : 1.0);
1310 
1311  h3.SetBinContent(ix, iy, iz, h3.GetBinContent(ix,iy,iz) / w);
1312 
1313  if (use_error) {
1314  h3.SetBinError(ix, iy, iz, h3.GetBinError(ix,iy,iz) / w);
1315  }
1316  }
1317  }
1318  }
1319  }
1320 
1321  } else if (Z) {
1322 
1323  for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1324  for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1325 
1326  if (normalise) {
1327 
1328  W = 0.0;
1329 
1330  for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1331  W += h3.GetBinContent(ix,iy,iz);
1332  }
1333  }
1334 
1335  if (W != 0.0) {
1336 
1337  for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1338 
1339  const Double_t w = W * (bin_width ? h3.GetZaxis()->GetBinWidth(iz) : 1.0);
1340 
1341  h3.SetBinContent(ix, iy, iz, h3.GetBinContent(ix,iy,iz) / w);
1342 
1343  if (use_error) {
1344  h3.SetBinError(ix, iy, iz, h3.GetBinError(ix,iy,iz) / w);
1345  }
1346  }
1347  }
1348  }
1349  }
1350  }
1351  }
1352 
1353 
1354  /**
1355  * Set limits of TGraph.
1356  *
1357  * \param g1 graph
1358  */
1359  inline void setLimits(TGraph& g1)
1360  {
1361  using namespace std;
1362 
1363  Double_t ymin = numeric_limits<Double_t>::max();
1364  Double_t ymax = numeric_limits<Double_t>::lowest();
1365 
1366  for (Int_t i = 0; i != g1.GetN(); ++i) {
1367 
1368  const Double_t y = g1.GetY()[i];
1369 
1370  if (y > ymax) { ymax = y; }
1371  if (y < ymin) { ymin = y; }
1372  }
1373 
1374  g1.SetMinimum(ymin);
1375  g1.SetMaximum(ymax);
1376  }
1377 
1378 
1379  /**
1380  * Set limits of TGraph2D.
1381  *
1382  * \param g2 graph
1383  */
1384  inline void setLimits(TGraph2D& g2)
1385  {
1386  using namespace std;
1387 
1388  Double_t zmin = numeric_limits<Double_t>::max();
1389  Double_t zmax = numeric_limits<Double_t>::lowest();
1390 
1391  for (Int_t i = 0; i != g2.GetN(); ++i) {
1392 
1393  const Double_t z = g2.GetZ()[i];
1394 
1395  if (z > zmax) { zmax = z; }
1396  if (z < zmin) { zmin = z; }
1397  }
1398 
1399  g2.SetMinimum(zmin);
1400  g2.SetMaximum(zmax);
1401  }
1402 
1403 
1404  /**
1405  * Set axis range.
1406  *
1407  * \param xmin lower limit (I/O)
1408  * \param xmax upper limit (I/O)
1409  * \param logx logarithmic
1410  */
1411  inline void setRange(double& xmin,
1412  double& xmax,
1413  const bool logx)
1414  {
1415  if (logx) {
1416  xmin = log(xmin);
1417  xmax = log(xmax);
1418  }
1419 
1420  double dx = (xmax - xmin) * 0.1;
1421 
1422  if (logx || xmin >= dx || xmin < 0.0)
1423  xmin -= dx;
1424  else
1425  xmin = 0.0;
1426 
1427  xmax += dx;
1428 
1429  if (logx) {
1430  xmin = exp(xmin);
1431  xmax = exp(xmax);
1432  }
1433  }
1434 
1435 
1436  /**
1437  * Set axis with PMT address labels.
1438  *
1439  * This methods sets the labels of the given axis to the sorted values of the PMT ring and position.\n
1440  * It should normally be called before filling of the corresponding histogram.\n
1441  * The filling should then be made with the textual representation of the PMT ring and position (i.e.\ JDETECTOR::JPMTPhysicalAddress::toString).
1442  *
1443  * Alternatively, the filling can be made with the index of the PMT in the address map of the corresponding module
1444  * (i.e.\ JDETECTOR::JModuleAddressMap::getIndex).\n
1445  * In that case, the method can be called before or after filling of the histogram.
1446  *
1447  * \param axis axis
1448  * \param memo module address map
1449  */
1450  inline void setAxisLabels(TAxis *axis, const JModuleAddressMap& memo)
1451  {
1452  using namespace JPP;
1453 
1454  if (axis->GetNbins() == (int) memo.size()) {
1455 
1456  for (int i = 0; i != axis->GetNbins(); ++i) {
1457 
1458  const JPMTPhysicalAddress& address = memo[i];
1459 
1460  axis->SetBinLabel(i + 1, address.toString().c_str());
1461  }
1462 
1463  } else {
1464 
1465  THROW(JValueOutOfRange, "Number of bins " << axis->GetNbins() << " != " << memo.size());
1466  }
1467  }
1468 
1469 
1470  /**
1471  * Set axis labels with PMT addresses.
1472  *
1473  * This methods sets the labels of the given axis to the sorted values of the PMT ring and position.\n
1474  * It should be called after filling of the corresponding histogram.\n
1475  * The filling should have been made with the PMT number (e.g.\ KM3NETDAQ::JDAQHit::getPMT).
1476  *
1477  * \param h1 histogram
1478  * \param axis axis <tt>(x, X, y, Y, z, Z)</tt>
1479  * \param memo module address map
1480  */
1481  inline void setAxisLabels(TH1& h1, const std::string axis, const JModuleAddressMap& memo)
1482  {
1483  using namespace JPP;
1484 
1485  TAxis* pax = NULL;
1486 
1487  if (axis == "x" || axis == "X")
1488  pax = h1.GetXaxis();
1489  else if (axis == "y" || axis == "Y")
1490  pax = h1.GetYaxis();
1491  else if (axis == "z" || axis == "Z")
1492  pax = h1.GetZaxis();
1493  else
1494  THROW(JValueOutOfRange, "Invalid axis " << axis);
1495 
1496  if (pax->GetNbins() == (int) memo.size()) {
1497 
1498  for (int i = 0; i != pax->GetNbins(); ++i) {
1499 
1500  const JPMTPhysicalAddress& address = memo.getAddressTranslator(i);
1501 
1502  pax->SetBinLabel(i + 1, address.toString().c_str());
1503  }
1504 
1505  h1.LabelsOption("a", axis.c_str()); // sort labels
1506 
1507  } else {
1508 
1509  THROW(JValueOutOfRange, "Number of bins " << pax->GetNbins() << " != " << memo.size());
1510  }
1511  }
1512 
1513 
1514  /**
1515  * Check if given key corresponds to a TObject.
1516  *
1517  * \param key ROOT key
1518  * \return true if given key corresponds to a TObject; else false
1519  */
1520  inline bool isTObject(const TKey* key)
1521  {
1522  return (key != NULL && TClass::GetClass(key->GetClassName())->IsTObject());
1523  }
1524 }
1525 
1526 #endif
static const char *const sqrt()
Set contents to signed difference between squares.
const double xmax
Definition: JQuadrature.cc:24
data_type w[N+1][M+1]
Definition: JPolint.hh:867
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.
Exception for a functional operation.
Definition: JException.hh:142
static const char *const SAME_AS_INPUT()
Set name of output histogram to name of input histogram.
then usage $script[< detector identifier >< run range >]< QA/QCfile > nExample script to produce data quality plots nWhen a detector identifier and run range are data are downloaded from the database nand subsequently stored in the given QA QC file
Definition: JDataQuality.sh:19
char text[TEXT_SIZE]
Definition: elog.cc:72
then set_variable DIR else fatal Wrong number of arguments fi for INPUT_FILE in ls rt $DIR stage * log
#define THROW(JException_t, A)
Marco for throwing exception with std::ostream compatible message.
Definition: JException.hh:712
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.
Auxiliary class to handle file name, ROOT directory and object name.
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 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 using log10()).
const JPolynome f1(1.0, 2.0, 3.0)
Function.
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.
void setLogarithmicX(TList *list)
Make x-axis of objects in list logarithmic (e.g. after using log10()).
Double_t getResult(const TString &text, TObject *object=NULL)
Get result of given textual formula.
do set_variable OUTPUT_DIRECTORY $WORKDIR T
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.
then awk string
T pow(const T &x, const double y)
Power .
Definition: JMath.hh:97
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 setLogarithmicY(TList *list)
Make y-axis of objects in list logarithmic (e.g. after using log10()).
const double xmin
Definition: JQuadrature.cc:23
void setRange(double &xmin, double &xmax, const bool logx)
Set axis range.
const JPMTAddressTranslator & getAddressTranslator(const int tdc) const
Get PMT address translator.
then usage $script< input file >[option[primary[working directory]]] nWhere option can be N
Definition: JMuonPostfit.sh:40
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 JAcoustics sh $DETECTOR_ID source JAcousticsToolkit sh CHECK_EXIT_CODE typeset A EMITTERS get_tripods $WORKDIR tripod txt EMITTERS get_transmitters $WORKDIR transmitter txt EMITTERS for EMITTER in
Definition: JCanberra.sh:48
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:162
Exception for parsing value.
Definition: JException.hh:196
bool isTObject(const TKey *key)
Check if given key corresponds to a TObject.
no fit printf nominal n $STRING awk v X
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:128
do set_variable MODULE getModule a $WORKDIR detector_a datx L $STRING JEditDetector a $WORKDIR detector_a datx M $MODULE setz o $WORKDIR detector_a datx JEditDetector a $WORKDIR detector_b datx M $MODULE setz o $WORKDIR detector_b datx done echo Output stored at $WORKDIR detector_a datx and $WORKDIR tripod_a txt JDrawDetector2D a $WORKDIR detector_a datx a $WORKDIR detector_b datx L BL o detector $FORMAT $BATCH JDrawDetector2D T $WORKDIR tripod_a txt T $WORKDIR tripod_b txt L BL o tripod $FORMAT $BATCH JCompareDetector a $WORKDIR detector_a datx b $WORKDIR detector_b datx o $WORKDIR abc root &dev null for KEY in X Y Z
Exception for accessing a value in a collection that is outside of its range.
Definition: JException.hh:178
TString getLogarithmic(const TString &formula, const char parameter)
Make given parameter in formula logarithmic (e.g. after using log10()).
static const char *const TIMESTAMP
Time stamp of earliest UTC time.
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"` source JAcousticsToolkit.sh typeset -A TRIPODS get_tripods $WORKDIR/tripod.txt TRIPODS XMEAN
Double_t g1(const Double_t x)
Function.
Definition: JQuantiles.cc:25