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