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