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