Jpp master_rocky-44-g75b7c4f75
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(f1->GetXmin(),
609 f1->GetXmax());
610
611 *f1 = fn;
612
613 f1->SetRange(TMath::Power(10.0,f1->GetXmin()),
614 TMath::Power(10.0,f1->GetXmax()));
615 }
616 }
617
618
619 /**
620 * Make y-axis of given function logarithmic (e.g.\ after using <tt>log10()</tt>).
621 *
622 * \param f1 function
623 */
624 inline void setLogarithmicY(TF1* f1)
625 {
626 if (f1 != NULL) {
627
628 TString buffer = f1->GetExpFormula();
629
630 buffer = "pow(10.0, " + buffer + ")";
631
632 TF1 fn(f1->GetName(), buffer);
633
634 copy(*f1, fn);
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 x-axis and associated functions of given histogram logarithmic (e.g.\ after using <tt>log10()</tt>).
715 *
716 * \param h2 histogram
717 */
718 inline void setLogarithmicX(TH2* h2)
719 {
720 using namespace std;
721
722 if (h2 != NULL) {
723
724 setLogarithmicX<TF2>(h2->GetListOfFunctions());
725
726 setLogarithmic(h2->GetXaxis());
727 }
728 }
729
730
731 /**
732 * Make y-axis and associated functions of given histogram logarithmic (e.g.\ after using <tt>log10()</tt>).
733 *
734 * \param h2 histogram
735 */
736 inline void setLogarithmicY(TH2* h2)
737 {
738 if (h2 != NULL) {
739
740 setLogarithmicY<TF2>(h2->GetListOfFunctions());
741
742 setLogarithmic(h2->GetYaxis());
743 }
744 }
745
746
747 /**
748 * Make x-axis and associated functions of given graph logarithmic (e.g.\ after using <tt>log10()</tt>).
749 *
750 * \param g1 graph
751 */
752 inline void setLogarithmicX(TGraph* g1)
753 {
754 if (g1 != NULL) {
755
756 setLogarithmicX<TF1>(g1->GetListOfFunctions());
757
758 for (Int_t i = 0; i != g1->GetN(); ++i) {
759 g1->GetX()[i] = pow(10.0, g1->GetX()[i]);
760 }
761 }
762 }
763
764
765 /**
766 * Make y-axis and associated functions of given graph logarithmic (e.g.\ after using <tt>log10()</tt>).
767 *
768 * \param g1 graph
769 */
770 inline void setLogarithmicY(TGraph* g1)
771 {
772 if (g1 != NULL) {
773
774 setLogarithmicY<TF1>(g1->GetListOfFunctions());
775
776 for (Int_t i = 0; i != g1->GetN(); ++i) {
777 g1->GetY()[i] = pow(10.0, g1->GetY()[i]);
778 }
779 }
780 }
781
782
783 /**
784 * Make x-axis and associated functions of given graph logarithmic (e.g.\ after using <tt>log10()</tt>).
785 *
786 * \param g1 graph
787 */
788 inline void setLogarithmicX(TGraphErrors* g1)
789 {
790 if (g1 != NULL) {
791
792 setLogarithmicX<TF1>(g1->GetListOfFunctions());
793
794 for (Int_t i = 0; i != g1->GetN(); ++i) {
795 g1->GetEX()[i] = pow(10.0, g1->GetX()[i] + g1->GetEX()[i]) - pow(10.0, g1->GetX()[i]);
796 g1->GetX() [i] = pow(10.0, g1->GetX()[i]);
797 }
798 }
799 }
800
801
802 /**
803 * Make y-axis and associated functions of given graph logarithmic (e.g.\ after using <tt>log10()</tt>).
804 *
805 * \param g1 graph
806 */
807 inline void setLogarithmicY(TGraphErrors* g1)
808 {
809 if (g1 != NULL) {
810
811 setLogarithmicY<TF1>(g1->GetListOfFunctions());
812
813 for (Int_t i = 0; i != g1->GetN(); ++i) {
814 g1->GetEY()[i] = pow(10.0, g1->GetY()[i] + g1->GetEY()[i]) - pow(10.0, g1->GetY()[i]);
815 g1->GetY() [i] = pow(10.0, g1->GetY()[i]);
816 }
817 }
818 }
819
820
821 /**
822 * Make x-axis and associated functions of given graph logarithmic (e.g.\ after using <tt>log10()</tt>).
823 *
824 * \param g2 graph
825 */
826 inline void setLogarithmicX(TGraph2D* g2)
827 {
828 if (g2 != NULL) {
829
830 setLogarithmicX<TF2>(g2->GetListOfFunctions());
831
832 for (Int_t i = 0; i != g2->GetN(); ++i) {
833 g2->GetX()[i] = pow(10.0, g2->GetX()[i]);
834 }
835 }
836 }
837
838
839 /**
840 * Make y-axis and associated functions of given graph logarithmic (e.g.\ after using <tt>log10()</tt>).
841 *
842 * \param g2 graph
843 */
844 inline void setLogarithmicY(TGraph2D* g2)
845 {
846 if (g2 != NULL) {
847
848 setLogarithmicY<TF2>(g2->GetListOfFunctions());
849
850 for (Int_t i = 0; i != g2->GetN(); ++i) {
851 g2->GetY()[i] = pow(10.0, g2->GetY()[i]);
852 }
853 }
854 }
855
856
857 /**
858 * Make x-axis and associated functions of given graph logarithmic (e.g.\ after using <tt>log10()</tt>).
859 *
860 * \param g2 graph
861 */
862 inline void setLogarithmicX(TGraph2DErrors* g2)
863 {
864 if (g2 != NULL) {
865
866 setLogarithmicX<TF2>(g2->GetListOfFunctions());
867
868 for (Int_t i = 0; i != g2->GetN(); ++i) {
869 g2->GetEX()[i] = pow(10.0, g2->GetX()[i] + g2->GetEX()[i]) - pow(10.0, g2->GetX()[i]);
870 g2->GetX() [i] = pow(10.0, g2->GetX()[i]);
871 }
872 }
873 }
874
875
876 /**
877 * Make y-axis and associated functions of given graph logarithmic (e.g.\ after using <tt>log10()</tt>).
878 *
879 * \param g2 graph
880 */
881 inline void setLogarithmicY(TGraph2DErrors* g2)
882 {
883 if (g2 != NULL) {
884
885 setLogarithmicY<TF2>(g2->GetListOfFunctions());
886
887 for (Int_t i = 0; i != g2->GetN(); ++i) {
888 g2->GetEY()[i] = pow(10.0, g2->GetY()[i] + g2->GetEY()[i]) - pow(10.0, g2->GetY()[i]);
889 g2->GetY() [i] = pow(10.0, g2->GetY()[i]);
890 }
891 }
892 }
893
894
895 /**
896 * Make x-axis of given multi-graph logarithmic (e.g.\ after using <tt>log10()</tt>).
897 *
898 * \param gn multi graph
899 */
900 inline void setLogarithmicX(TMultiGraph* gn)
901 {
902 if (gn != NULL) {
903 setLogarithmicX<TGraph>(gn->GetListOfGraphs());
904 }
905 }
906
907
908 /**
909 * Make y-axis of given multi-graph logarithmic (e.g.\ after using <tt>log10()</tt>).
910 *
911 * \param gn multi graph
912 */
913 inline void setLogarithmicY(TMultiGraph* gn)
914 {
915 if (gn != NULL) {
916 setLogarithmicY<TGraph>(gn->GetListOfGraphs());
917 }
918 }
919
920
921 /**
922 * Make x-axis of given line logarithmic (e.g.\ after using <tt>log10()</tt>).
923 *
924 * \param line line
925 */
926 inline void setLogarithmicX(TLine* line)
927 {
928 if (line != NULL) {
929 line->SetX1(pow(10.0, line->GetX1()));
930 line->SetX2(pow(10.0, line->GetX2()));
931 }
932 }
933
934
935 /**
936 * Make y-axis of given line logarithmic (e.g.\ after using <tt>log10()</tt>).
937 *
938 * \param line line
939 */
940 inline void setLogarithmicY(TLine* line)
941 {
942 if (line != NULL) {
943 line->SetY1(pow(10.0, line->GetY1()));
944 line->SetY2(pow(10.0, line->GetY2()));
945 }
946 }
947
948
949 /**
950 * Make x-axis of given ellipse logarithmic (e.g.\ after using <tt>log10()</tt>).
951 *
952 * \param ellipse ellipse
953 */
954 inline void setLogarithmicX(TEllipse* ellipse)
955 {
956 THROW(JFunctionalException, "Operation setLogarithmicX on TEllipse not allowed.");
957 }
958
959
960 /**
961 * Make y-axis of given ellipse logarithmic (e.g.\ after using <tt>log10()</tt>).
962 *
963 * \param ellipse ellipse
964 */
965 inline void setLogarithmicY(TEllipse* ellipse)
966 {
967 THROW(JFunctionalException, "Operation setLogarithmicY on TEllipse not allowed.");
968 }
969
970
971 /**
972 * Make x-axis of objects in list logarithmic (e.g.\ after using <tt>log10()</tt>).
973 *
974 * \param list list
975 */
976 template<class T>
977 inline void setLogarithmicX(TList* list)
978 {
979 for (TIter i(list); T* p = dynamic_cast<T*>(i.Next()); ) {
981 }
982 }
983
984
985 /**
986 * Make y-axis of objects in list logarithmic (e.g.\ after using <tt>log10()</tt>).
987 *
988 * \param list list
989 */
990 template<class T>
991 inline void setLogarithmicY(TList* list)
992 {
993 for (TIter i(list); T* p = dynamic_cast<T*>(i.Next()); ) {
995 }
996 }
997
998
999 /**
1000 * Convert 1D histogram to PDF.
1001 *
1002 * Possible options are:
1003 * - N normalise to histogram contents;
1004 * - W divide by bin width;
1005 * - E convert also bin errors.
1006 *
1007 * \param h1 histogram
1008 * \param option option
1009 * \param factor scaling factor
1010 */
1011 inline void convertToPDF(TH1& h1, const std::string& option = "NW", const double factor = 1.0)
1012 {
1013 using namespace std;
1014
1015 const bool normalise = (option.find('N') != string::npos || option.find('n') != string::npos);
1016 const bool bin_width = (option.find('W') != string::npos || option.find('w') != string::npos);
1017 const bool use_error = (option.find('E') != string::npos || option.find('e') != string::npos);
1018
1019 Double_t W = 1.0;
1020
1021 if (normalise) {
1022
1023 W = 0.0;
1024
1025 for (Int_t i = 1; i <= h1.GetXaxis()->GetNbins(); ++i) {
1026 W += h1.GetBinContent(i);
1027 }
1028 }
1029
1030 if (W != 0.0) {
1031
1032 for (Int_t i = 1; i <= h1.GetXaxis()->GetNbins(); ++i) {
1033
1034 const Double_t w = W * (bin_width ? h1.GetXaxis()->GetBinWidth(i) : 1.0);
1035
1036 h1.SetBinContent(i, h1.GetBinContent(i) * factor / w);
1037
1038 if (use_error) {
1039 h1.SetBinError(i, h1.GetBinError(i) * factor / w);
1040 }
1041 }
1042 }
1043 }
1044
1045
1046 /**
1047 * Convert 2D histogram to PDF.
1048 *
1049 * Possible options are:
1050 * - N normalise to histogram contents;
1051 * - X convert x-axis to PDF;
1052 * - Y convert y-axis to PDF;
1053 * - W divide by bin width;
1054 * - E convert also bin errors.
1055 *
1056 * \param h2 histogram
1057 * \param option option
1058 * \param factor scaling factor
1059 */
1060 inline void convertToPDF(TH2& h2, const std::string& option = "NXYW", const double factor = 1.0)
1061 {
1062 using namespace std;
1063
1064 const bool normalise = (option.find('N') != string::npos || option.find('n') != string::npos);
1065 const bool X = (option.find('X') != string::npos || option.find('x') != string::npos);
1066 const bool Y = (option.find('Y') != string::npos || option.find('y') != string::npos);
1067 const bool bin_width = (option.find('W') != string::npos || option.find('w') != string::npos);
1068 const bool use_error = (option.find('E') != string::npos || option.find('e') != string::npos);
1069
1070 Double_t W = 1.0;
1071
1072 if (X && Y) {
1073
1074 if (normalise) {
1075
1076 W = 0.0;
1077
1078 for (Int_t ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
1079 for (Int_t iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
1080 W += h2.GetBinContent(ix,iy);
1081 }
1082 }
1083 }
1084
1085 if (W != 0.0) {
1086
1087 for (Int_t ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
1088 for (Int_t iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
1089
1090 const Double_t w = W * (bin_width ? h2.GetXaxis()->GetBinWidth(ix) * h2.GetYaxis()->GetBinWidth(iy) : 1.0);
1091
1092 h2.SetBinContent(ix, iy, h2.GetBinContent(ix,iy) * factor / w);
1093
1094 if (use_error) {
1095 h2.SetBinError(ix, iy, h2.GetBinError(ix,iy) * factor / w);
1096 }
1097 }
1098 }
1099 }
1100
1101 } else if (X) {
1102
1103 for (Int_t iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
1104
1105 if (normalise) {
1106
1107 W = 0.0;
1108
1109 for (Int_t ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
1110 W += h2.GetBinContent(ix,iy);
1111 }
1112 }
1113
1114 if (W != 0.0) {
1115
1116 for (Int_t ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
1117
1118 const Double_t w = W * (bin_width ? h2.GetXaxis()->GetBinWidth(ix) : 1.0);
1119
1120 h2.SetBinContent(ix, iy, h2.GetBinContent(ix,iy) * factor / w);
1121
1122 if (use_error) {
1123 h2.SetBinError(ix, iy, h2.GetBinError(ix,iy) * factor / w);
1124 }
1125 }
1126 }
1127 }
1128
1129 } else if (Y) {
1130
1131 for (Int_t ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
1132
1133 if (normalise) {
1134
1135 W = 0.0;
1136
1137 for (Int_t iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
1138 W += h2.GetBinContent(ix,iy);
1139 }
1140 }
1141
1142 if (W != 0.0) {
1143
1144 for (Int_t iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
1145
1146 const Double_t w = W * (bin_width ? h2.GetYaxis()->GetBinWidth(iy) : 1.0);
1147
1148 h2.SetBinContent(ix, iy, h2.GetBinContent(ix,iy) / w);
1149
1150 if (use_error) {
1151 h2.SetBinError(ix, iy, h2.GetBinError(ix,iy) / w);
1152 }
1153 }
1154 }
1155 }
1156 }
1157 }
1158
1159
1160 /**
1161 * Convert 3D histogram to PDF.
1162 *
1163 * Possible options are:
1164 * - N normalise to histogram contents;
1165 * - X convert x-axis to PDF;
1166 * - Y convert y-axis to PDF;
1167 * - Z convert z-axis to PDF;
1168 * - W divide by bin width;
1169 * - E convert also bin errors.
1170 *
1171 * \param h3 histogram
1172 * \param option option
1173 * \param factor scaling factor
1174 */
1175 inline void convertToPDF(TH3& h3, const std::string& option = "NXYW", const double factor = 1.0)
1176 {
1177 using namespace std;
1178
1179 const bool normalise = (option.find('N') != string::npos || option.find('n') != string::npos);
1180 const bool X = (option.find('X') != string::npos || option.find('x') != string::npos);
1181 const bool Y = (option.find('Y') != string::npos || option.find('y') != string::npos);
1182 const bool Z = (option.find('Z') != string::npos || option.find('z') != string::npos);
1183 const bool bin_width = (option.find('W') != string::npos || option.find('w') != string::npos);
1184 const bool use_error = (option.find('E') != string::npos || option.find('e') != string::npos);
1185
1186 Double_t W = 1.0;
1187
1188 if (X && Y && Z) {
1189
1190 if (normalise) {
1191
1192 W = 0.0;
1193
1194 for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1195 for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1196 for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1197 W += h3.GetBinContent(ix,iy,iz);
1198 }
1199 }
1200 }
1201 }
1202
1203 if (W != 0.0) {
1204
1205 for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1206 for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1207 for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1208
1209 const Double_t w = W * (bin_width ? h3.GetXaxis()->GetBinWidth(ix) * h3.GetYaxis()->GetBinWidth(iy) * h3.GetZaxis()->GetBinWidth(iz) : 1.0);
1210
1211 h3.SetBinContent(ix, iy, iz, h3.GetBinContent(ix,iy,iz) * factor / w);
1212
1213 if (use_error) {
1214 h3.SetBinError(ix, iy, iz, h3.GetBinError(ix,iy,iz) * factor / w);
1215 }
1216 }
1217 }
1218 }
1219 }
1220
1221 } else if (X && Z) {
1222
1223 for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1224
1225 if (normalise) {
1226
1227 W = 0.0;
1228
1229 for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1230 for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1231 W += h3.GetBinContent(ix,iy,iz);
1232 }
1233 }
1234 }
1235
1236 if (W != 0.0) {
1237
1238 for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1239 for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1240
1241 const Double_t w = W * (bin_width ? h3.GetXaxis()->GetBinWidth(ix) * h3.GetZaxis()->GetBinWidth(iz) : 1.0);
1242
1243 h3.SetBinContent(ix, iy, iz, h3.GetBinContent(ix,iy,iz) * factor / w);
1244
1245 if (use_error) {
1246 h3.SetBinError(ix, iy, iz, h3.GetBinError(ix,iy,iz) * factor / w);
1247 }
1248 }
1249 }
1250 }
1251 }
1252
1253 } else if (Y && Z) {
1254
1255 for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1256
1257 if (normalise) {
1258
1259 W = 0.0;
1260
1261 for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1262 for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1263 W += h3.GetBinContent(ix,iy,iz);
1264 }
1265 }
1266 }
1267
1268 if (W != 0.0) {
1269
1270 for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1271 for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1272
1273 const Double_t w = W * (bin_width ? h3.GetYaxis()->GetBinWidth(iy) * h3.GetZaxis()->GetBinWidth(iz) : 1.0);
1274
1275 h3.SetBinContent(ix, iy, iz, h3.GetBinContent(ix,iy,iz) * factor / w);
1276
1277 if (use_error) {
1278 h3.SetBinError(ix, iy, iz, h3.GetBinError(ix,iy,iz) * factor / w);
1279 }
1280 }
1281 }
1282 }
1283 }
1284
1285 } else if (X && Y) {
1286
1287 for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1288
1289 if (normalise) {
1290
1291 W = 0.0;
1292
1293 for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1294 for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1295 W += h3.GetBinContent(ix,iy,iz);
1296 }
1297 }
1298 }
1299
1300 if (W != 0.0) {
1301
1302 for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1303 for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1304
1305 const Double_t w = W * (bin_width ? h3.GetXaxis()->GetBinWidth(ix) * h3.GetYaxis()->GetBinWidth(iy) : 1.0);
1306
1307 h3.SetBinContent(ix, iy, iz, h3.GetBinContent(ix,iy,iz) * factor / w);
1308
1309 if (use_error) {
1310 h3.SetBinError(ix, iy, iz, h3.GetBinError(ix,iy,iz) * factor / w);
1311 }
1312 }
1313 }
1314 }
1315 }
1316
1317 } else if (X) {
1318
1319 for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1320 for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1321
1322 if (normalise) {
1323
1324 W = 0.0;
1325
1326 for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1327 W += h3.GetBinContent(ix,iy,iz);
1328 }
1329 }
1330
1331 if (W != 0.0) {
1332
1333 for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1334
1335 const Double_t w = W * (bin_width ? h3.GetXaxis()->GetBinWidth(ix) : 1.0);
1336
1337 h3.SetBinContent(ix, iy, iz, h3.GetBinContent(ix,iy,iz) * factor / w);
1338
1339 if (use_error) {
1340 h3.SetBinError(ix, iy, iz, h3.GetBinError(ix,iy,iz) * factor / w);
1341 }
1342 }
1343 }
1344 }
1345 }
1346
1347 } else if (Y) {
1348
1349 for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1350 for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1351
1352 if (normalise) {
1353
1354 W = 0.0;
1355
1356 for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1357 W += h3.GetBinContent(ix,iy,iz);
1358 }
1359 }
1360
1361 if (W != 0.0) {
1362
1363 for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1364
1365 const Double_t w = W * (bin_width ? h3.GetYaxis()->GetBinWidth(iy) : 1.0);
1366
1367 h3.SetBinContent(ix, iy, iz, h3.GetBinContent(ix,iy,iz) / w);
1368
1369 if (use_error) {
1370 h3.SetBinError(ix, iy, iz, h3.GetBinError(ix,iy,iz) / w);
1371 }
1372 }
1373 }
1374 }
1375 }
1376
1377 } else if (Z) {
1378
1379 for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1380 for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1381
1382 if (normalise) {
1383
1384 W = 0.0;
1385
1386 for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1387 W += h3.GetBinContent(ix,iy,iz);
1388 }
1389 }
1390
1391 if (W != 0.0) {
1392
1393 for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1394
1395 const Double_t w = W * (bin_width ? h3.GetZaxis()->GetBinWidth(iz) : 1.0);
1396
1397 h3.SetBinContent(ix, iy, iz, h3.GetBinContent(ix,iy,iz) / w);
1398
1399 if (use_error) {
1400 h3.SetBinError(ix, iy, iz, h3.GetBinError(ix,iy,iz) / w);
1401 }
1402 }
1403 }
1404 }
1405 }
1406 }
1407 }
1408
1409
1410 /**
1411 * Set limits of TGraph.
1412 *
1413 * \param g1 graph
1414 */
1415 inline void setLimits(TGraph& g1)
1416 {
1417 using namespace std;
1418
1419 Double_t ymin = numeric_limits<Double_t>::max();
1420 Double_t ymax = numeric_limits<Double_t>::lowest();
1421
1422 for (Int_t i = 0; i != g1.GetN(); ++i) {
1423
1424 const Double_t y = g1.GetY()[i];
1425
1426 if (y > ymax) { ymax = y; }
1427 if (y < ymin) { ymin = y; }
1428 }
1429
1430 g1.SetMinimum(ymin);
1431 g1.SetMaximum(ymax);
1432 }
1433
1434
1435 /**
1436 * Set limits of TGraph2D.
1437 *
1438 * \param g2 graph
1439 */
1440 inline void setLimits(TGraph2D& g2)
1441 {
1442 using namespace std;
1443
1444 Double_t zmin = numeric_limits<Double_t>::max();
1445 Double_t zmax = numeric_limits<Double_t>::lowest();
1446
1447 for (Int_t i = 0; i != g2.GetN(); ++i) {
1448
1449 const Double_t z = g2.GetZ()[i];
1450
1451 if (z > zmax) { zmax = z; }
1452 if (z < zmin) { zmin = z; }
1453 }
1454
1455 g2.SetMinimum(zmin);
1456 g2.SetMaximum(zmax);
1457 }
1458
1459
1460 /**
1461 * Set axis range.
1462 *
1463 * \param xmin lower limit (I/O)
1464 * \param xmax upper limit (I/O)
1465 * \param logx logarithmic
1466 */
1467 inline void setRange(double& xmin,
1468 double& xmax,
1469 const bool logx)
1470 {
1471 if (logx) {
1472 xmin = log(xmin);
1473 xmax = log(xmax);
1474 }
1475
1476 double dx = (xmax - xmin) * 0.1;
1477
1478 if (logx || xmin >= dx || xmin < 0.0)
1479 xmin -= dx;
1480 else
1481 xmin = 0.0;
1482
1483 xmax += dx;
1484
1485 if (logx) {
1486 xmin = exp(xmin);
1487 xmax = exp(xmax);
1488 }
1489 }
1490
1491
1492 /**
1493 * Set axis with PMT address labels.
1494 *
1495 * This methods sets the labels of the given axis to the sorted values of the PMT ring and position.\n
1496 * It should normally be called before filling of the corresponding histogram.\n
1497 * The filling should then be made with the textual representation of the PMT ring and position (i.e.\ JDETECTOR::JPMTPhysicalAddress::toString).
1498 *
1499 * Alternatively, the filling can be made with the index of the PMT in the address map of the corresponding module
1500 * (i.e.\ JDETECTOR::JModuleAddressMap::getIndex).\n
1501 * In that case, the method can be called before or after filling of the histogram.
1502 *
1503 * \param axis axis
1504 * \param memo module address map
1505 */
1506 inline void setAxisLabels(TAxis *axis, const JModuleAddressMap& memo)
1507 {
1508 using namespace JPP;
1509
1510 if (axis->GetNbins() == (int) memo.size()) {
1511
1512 for (int i = 0; i != axis->GetNbins(); ++i) {
1513
1514 const JPMTPhysicalAddress& address = memo[i];
1515
1516 axis->SetBinLabel(i + 1, address.toString().c_str());
1517 }
1518
1519 } else {
1520
1521 THROW(JValueOutOfRange, "Number of bins " << axis->GetNbins() << " != " << memo.size());
1522 }
1523 }
1524
1525
1526 /**
1527 * Set axis labels with PMT addresses.
1528 *
1529 * This methods sets the labels of the given axis to the sorted values of the PMT ring and position.\n
1530 * It should be called after filling of the corresponding histogram.\n
1531 * The filling should have been made with the PMT number (e.g.\ KM3NETDAQ::JDAQHit::getPMT).
1532 *
1533 * \param h1 histogram
1534 * \param axis axis <tt>(x, X, y, Y, z, Z)</tt>
1535 * \param memo module address map
1536 */
1537 inline void setAxisLabels(TH1& h1, const std::string axis, const JModuleAddressMap& memo)
1538 {
1539 using namespace JPP;
1540
1541 TAxis* pax = NULL;
1542
1543 if (axis == "x" || axis == "X")
1544 pax = h1.GetXaxis();
1545 else if (axis == "y" || axis == "Y")
1546 pax = h1.GetYaxis();
1547 else if (axis == "z" || axis == "Z")
1548 pax = h1.GetZaxis();
1549 else
1550 THROW(JValueOutOfRange, "Invalid axis " << axis);
1551
1552 if (pax->GetNbins() == (int) memo.size()) {
1553
1554 for (int i = 0; i != pax->GetNbins(); ++i) {
1555
1556 const JPMTPhysicalAddress& address = memo.getAddressTranslator(i);
1557
1558 pax->SetBinLabel(i + 1, address.toString().c_str());
1559 }
1560
1561 h1.LabelsOption("a", axis.c_str()); // sort labels
1562
1563 } else {
1564
1565 THROW(JValueOutOfRange, "Number of bins " << pax->GetNbins() << " != " << memo.size());
1566 }
1567 }
1568
1569
1570 /**
1571 * Check if given key corresponds to a TObject.
1572 *
1573 * \param key ROOT key
1574 * \return true if given key corresponds to a TObject; else false
1575 */
1576 inline bool isTObject(const TKey* key)
1577 {
1578 return (key != NULL && TClass::GetClass(key->GetClassName())->IsTObject());
1579 }
1580
1581
1582 /**
1583 * Get first object of which name matches given reguar expression.
1584 *
1585 * \param ls pointer to list of objects
1586 * \param name regular expression
1587 * \return pointer to object (maybe NULL)
1588 */
1589 inline TObject* getObject(TList* ls, const char* const name)
1590 {
1591 const TRegexp regexp(name);
1592
1593 TIter iter(ls);
1594
1595 for (TObject* p1; (p1 = (TObject*) iter.Next()) != NULL; ) {
1596
1597 const TString tag(p1->GetName());
1598
1599 if (tag.Contains(regexp)) {
1600 return p1;
1601 }
1602 }
1603
1604 return NULL;
1605 }
1606
1607
1608 /**
1609 * Get function.
1610 *
1611 * \param h1 histogram
1612 * \param fcn function name
1613 * \return pointer to function (maybe NULL)
1614 */
1615 inline TF1* getFunction(TH1* h1, const char* const fcn)
1616 {
1617 return (TF1*) getObject(h1->GetListOfFunctions(), fcn);
1618 }
1619
1620
1621 /**
1622 * Get function.
1623 *
1624 * \param g1 graph
1625 * \param fcn function name
1626 * \return pointer to function (maybe NULL)
1627 */
1628 inline TF1* getFunction(TGraph* g1, const char* const fcn)
1629 {
1630 return (TF1*) getObject(g1->GetListOfFunctions(), fcn);
1631 }
1632
1633
1634 /**
1635 * Get function.
1636 *
1637 * \param g2 graph
1638 * \param fcn function name
1639 * \return pointer to function (maybe NULL)
1640 */
1641 inline TF1* getFunction(TGraph2D* g2, const char* const fcn)
1642 {
1643 return (TF1*) getObject(g2->GetListOfFunctions(), fcn);
1644 }
1645}
1646
1647#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.