Jpp  18.5.2
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JRootToolkit.hh
Go to the documentation of this file.
1 #ifndef __JROOT__JROOTTOOLKIT__
2 #define __JROOT__JROOTTOOLKIT__
3 
4 #include <string>
5 #include <istream>
6 #include <ostream>
7 #include <limits>
8 #include <cctype>
9 #include <vector>
10 
11 #pragma GCC diagnostic push
12 #pragma GCC diagnostic ignored "-Wall"
13 #include "TString.h"
14 #include "TObjString.h"
15 #include "TRegexp.h"
16 #include "TPRegexp.h"
17 #include "TFormula.h"
18 #include "TFile.h"
19 #include "TStreamerInfo.h"
20 #include "TList.h"
21 #include "TIterator.h"
22 #include "TF1.h"
23 #include "TH1.h"
24 #include "TH2.h"
25 #include "TGraph.h"
26 #include "TGraphErrors.h"
27 #include "TGraph2D.h"
28 #include "TGraph2DErrors.h"
29 #include "TMultiGraph.h"
30 #include "TNtuple.h"
31 #pragma GCC diagnostic pop
32 
33 #include "JLang/JType.hh"
34 #include "JLang/JLangToolkit.hh"
35 
36 #include "Jeep/JPrint.hh"
37 
38 #include "JROOT/JRootFile.hh"
39 #include "JROOT/JRootDictionary.hh"
40 #include "JROOT/JRootPrinter.hh"
41 
42 /**
43  * \author mdejong, mlincett
44  */
45 
46 namespace JROOT {}
47 namespace JPP { using namespace JROOT; }
48 
49 namespace JROOT {
50 
51  using JLANG::JType;
52 
53 
54  /**
55  * Get ROOT name of given data type.
56  *
57  * \return name of object to be read
58  */
59  template<class T>
60  inline const char* getName()
61  {
62  return getName(JType<T>());
63  }
64 
65 
66  /**
67  * Get ROOT name of given data type.
68  *
69  * \param type data type
70  * \return name of object to be read
71  */
72  template<class T>
73  inline const char* getName(const JType<T>& type)
74  {
75  return T::Class_Name();
76  }
77 
78 
79  /**
80  * Detach TH1 object and optionally reset contents.
81  *
82  * \param object ROOT TH1 object
83  * \param reset reset contents if true
84  * \return true if successful; else false
85  */
86  inline bool resetObject(TH1* object, const bool reset = false)
87  {
88  if (object != NULL) {
89 
90  object->SetDirectory(0);
91 
92  if (reset) {
93  object->Reset();
94  }
95 
96  return true;
97  }
98 
99  return false;
100  }
101 
102 
103  /**
104  * Detach TGraph object and optionally reset contents.
105  *
106  * \param object ROOT TGraph object
107  * \param reset reset contents if true
108  * \return true if successful; else false
109  */
110  inline bool resetObject(TGraph* object, const bool reset = false)
111  {
112  if (object != NULL) {
113 
114  if (reset) {
115  for (int i = 0; i != object->GetN(); ++i) {
116  object->SetPoint(i, 0.0, 0.0);
117  }
118  }
119 
120  return true;
121  }
122 
123  return false;
124  }
125 
126 
127  /**
128  * Detach TGraphErrors object and optionally reset contents.
129  *
130  * \param object ROOT TGraphErrors object
131  * \param reset reset contents if true
132  * \return true if successful; else false
133  */
134  inline bool resetObject(TGraphErrors* object, const bool reset = false)
135  {
136  if (object != NULL) {
137 
138  if (reset) {
139 
140  for (int i = 0; i != object->GetN(); ++i) {
141  object->SetPoint (i, 0.0, 0.0);
142  object->SetPointError(i, 0.0, 0.0);
143  }
144  }
145 
146  return true;
147  }
148 
149  return false;
150  }
151 
152 
153  /**
154  * Detach TGraph2D object and optionally reset contents.
155  *
156  * \param object ROOT TGraph2D object
157  * \param reset reset contents if true
158  * \return true if successful; else false
159  */
160  inline bool resetObject(TGraph2D* object, const bool reset = false)
161  {
162  if (object != NULL) {
163 
164  if (reset) {
165  for (int i = 0; i != object->GetN(); ++i) {
166  object->SetPoint(i, 0.0, 0.0, 0.0);
167  }
168  }
169 
170  return true;
171  }
172 
173  return false;
174  }
175 
176 
177  /**
178  * Detach TGraph2DErrors object and optionally reset contents.
179  *
180  * \param object ROOT TGraph2DErrors object
181  * \param reset reset contents if true
182  * \return true if successful; else false
183  */
184  inline bool resetObject(TGraph2DErrors* object, const bool reset = false)
185  {
186  if (object != NULL) {
187 
188  if (reset) {
189  for (int i = 0; i != object->GetN(); ++i) {
190  object->SetPoint (i, 0.0, 0.0, 0.0);
191  object->SetPointError(i, 0.0, 0.0, 0.0);
192  }
193  }
194 
195  return true;
196  }
197 
198  return false;
199  }
200 
201 
202  /**
203  * Detach TMultiGraph object and optionally reset contents.
204  *
205  * \param object ROOT TMultiGraph object
206  * \param reset reset contents if true
207  * \return true if successful; else false
208  */
209  inline bool resetObject(TMultiGraph* object, const bool reset = false)
210  {
211  if (object != NULL) {
212 
213  bool status = true;
214 
215  if (reset) {
216 
217  TIter next(object->GetListOfGraphs());
218 
219  while (TGraph* graph = (TGraph*)next()) {
220 
221  if (!(resetObject(dynamic_cast<TGraph2DErrors*>(graph), reset) ||
222  resetObject(dynamic_cast<TGraph2D*> (graph), reset)) &&
223  !(resetObject(dynamic_cast<TGraphErrors*> (graph), reset) ||
224  resetObject (graph, reset)) ) {
225 
226  status = false;
227  }
228  }
229  }
230 
231  return status;
232  }
233 
234  return false;
235  }
236 
237 
238  /**
239  * Detach TTree object and optionally reset contents.
240  *
241  * \param object ROOT TTree object
242  * \param reset reset contents if true
243  * \return true if successful; else false
244  */
245  inline bool resetObject(TTree* object, const bool reset = false)
246  {
247  if (object != NULL) {
248 
249  object->SetDirectory(0);
250 
251  if (reset) {
252  object->Reset();
253  }
254 
255  return true;
256  }
257 
258  return false;
259  }
260 
261 
262  /**
263  * Add point to TGraph.
264  *
265  * \param g1 pointer to valid ROOT TGraph object
266  * \param x x value
267  * \param y y value
268  */
269  inline void AddPoint(TGraph* g1,
270  const Double_t x,
271  const Double_t y)
272  {
273  const Int_t n = g1->GetN();
274 
275  g1->Set(n + 1);
276  g1->SetPoint(n, x, y);
277  }
278 
279 
280  /**
281  * Add point to TGraph2D.
282  *
283  * \param g1 pointer to valid ROOT TGraph object
284  * \param x x value
285  * \param y y value
286  * \param z z value
287  */
288  inline void AddPoint(TGraph2D* g1,
289  const Double_t x,
290  const Double_t y,
291  const Double_t z)
292  {
293  const Int_t n = g1->GetN();
294 
295  g1->Set(n + 1);
296  g1->SetPoint(n, x, y, z);
297  }
298 
299 
300  /**
301  * Add point to TGraphErrors.
302  *
303  * \param g1 pointer to valid ROOT TGraph object
304  * \param x x value
305  * \param y y value
306  * \param ex x error
307  * \param ey y error
308  */
309  inline void AddPoint(TGraphErrors* g1,
310  const Double_t x,
311  const Double_t y,
312  const Double_t ex,
313  const Double_t ey)
314  {
315  const Int_t n = g1->GetN();
316 
317  g1->Set(n + 1);
318  g1->SetPoint(n, x, y);
319  g1->SetPointError(n, ex, ey);
320  }
321 
322 
323  /**
324  * Add point to TGraph2DErrors.
325  *
326  * \param g1 pointer to valid ROOT TGraph object
327  * \param x x value
328  * \param y y value
329  * \param z z value
330  * \param ex x error
331  * \param ey y error
332  * \param ez z error
333  */
334  inline void AddPoint(TGraph2DErrors* g1,
335  const Double_t x,
336  const Double_t y,
337  const Double_t z,
338  const Double_t ex,
339  const Double_t ey,
340  const Double_t ez)
341  {
342  const Int_t n = g1->GetN();
343 
344  g1->Set(n + 1);
345  g1->SetPoint(n, x, y, z);
346  g1->SetPointError(n, ex, ey, ez);
347  }
348 
349 
350  /**
351  * Write object to ROOT file.
352  *
353  * \param file ROOT file
354  * \param object ROOT object
355  * \return ROOT file
356  */
357  inline TFile& operator<<(TFile& file, const TObject& object)
358  {
359  file.WriteTObject(&object);
360 
361  return file;
362  }
363 
364 
365  /**
366  * Get ROOT streamer information of class with given name.
367  * Note that the class name should include the name space, if any.
368  *
369  * \param file pointer to ROOT file
370  * \param name class name
371  * \return pointer to TStreamerInfo (NULL in case of error)
372  */
373  inline const TStreamerInfo* getStreamerInfo(TFile* file, const char* const name)
374  {
375  if (file != NULL && file->IsOpen()) {
376 
377  const TString buffer(name);
378 
379  TIter iter(file->GetStreamerInfoList());
380 
381  for (const TStreamerInfo* pStreamerInfo; (pStreamerInfo = (TStreamerInfo*) iter.Next()) != NULL; ) {
382  if (buffer == TString(pStreamerInfo->GetName())) {
383  return pStreamerInfo;
384  }
385  }
386  }
387 
388  return NULL;
389  }
390 
391 
392  /**
393  * Get ROOT streamer information of class with given name.
394  * Note that the class name should include the name space, if any.
395  *
396  * \param file_name file name
397  * \param name class name
398  * \return pointer to TStreamerInfo (NULL in case of error)
399  */
400  inline const TStreamerInfo* getStreamerInfo(const char* const file_name, const char* const name)
401  {
402  JRootInputFile file(file_name);
403 
404  return getStreamerInfo(file.getFile(), name);
405  }
406 
407 
408  /**
409  * Get ROOT streamer version of class with given name.
410  * Note that the class name should include the name space, if any.
411  *
412  * \param file pointer to ROOT file
413  * \param name class name
414  * \return streamer version; (-1 in case of error)
415  */
416  inline int getStreamerVersion(TFile* file, const char* const name)
417  {
418  const TStreamerInfo* pStreamerInfo = getStreamerInfo(file, name);
419 
420  if (pStreamerInfo != NULL)
421  return pStreamerInfo->GetClassVersion();
422  else
423  return -1;
424  }
425 
426 
427  /**
428  * Get ROOT streamer version of class with given name.
429  * Note that the class name should include the name space, if any.
430  *
431  * \param file_name file name
432  * \param name class name
433  * \return streamer version; (-1 in case of error)
434  */
435  inline int getStreamerVersion(const char* const file_name, const char* const name)
436  {
437  JRootInputFile file(file_name);
438 
439  return getStreamerVersion(file.getFile(), name);
440  }
441 
442 
443  /**
444  * Match a regular expression with given string and return the specified matched parts.
445  *
446  * If no matches are found corresponding to the specified index, the original string is returned.
447  *
448  * \param regexp regular expression
449  * \param string input string
450  * \param index index of matched parts (starting at 1)
451  * \return matched part of string
452  */
453  inline TString parse(const TPRegexp& regexp, const TString& string, const int index = 1)
454  {
455  TPRegexp buffer = regexp;
456  TObjArray* array = buffer.MatchS(string);
457  TString result = string;
458 
459  if (index - 1 < array->GetLast()) {
460  result = ((TObjString*) array->At(index))->GetName();
461  }
462 
463  delete array;
464 
465  return result;
466  }
467 
468 
469  /**
470  * Auxiliary data structure for a parameter index and its value.
471  */
473  /**
474  * Default constructor.
475  */
477  index(0),
478  value(0.0)
479  {}
480 
481 
482  /**
483  * Constructor.
484  *
485  * \param index parameter index
486  * \param value parameter value
487  */
488  JFitParameter_t(const Int_t index,
489  const Double_t value) :
490  index(index),
491  value(value)
492  {}
493 
494 
495  /**
496  * Type conversion operator
497  *
498  *
499  * \return index
500  */
501  inline operator Int_t() const
502  {
503  return index;
504  }
505 
506 
507  Int_t index;
508  Double_t value;
509  };
510 
511 
512  /**
513  * Set fit parameter.
514  *
515  * \param f1 fit function
516  * \param parameter parameter index and value
517  */
518  inline bool setParameter(TF1& f1, const JFitParameter_t& parameter)
519  {
520  if (parameter.index >= 0 && parameter.index < f1.GetNpar()) {
521 
522  f1.SetParameter(parameter.index, parameter.value);
523 
524  return true;
525 
526  } else {
527 
528  return false;
529  }
530  }
531 
532 
533  /**
534  * Fix fit parameter.
535  *
536  * \param f1 fit function
537  * \param parameter parameter index and value
538  */
539  inline bool fixParameter(TF1& f1, const JFitParameter_t& parameter)
540  {
541  if (parameter.index >= 0 && parameter.index < f1.GetNpar()) {
542 
543  f1.FixParameter(parameter.index, parameter.value);
544 
545  return true;
546 
547  } else {
548 
549  return false;
550  }
551  }
552 
553 
554  /**
555  * Release fit parameter.
556  *
557  * \param f1 fit function
558  * \param index parameter index
559  */
560  inline bool releaseParameter(TF1& f1, const Int_t index)
561  {
562  if (index >= 0 && index < f1.GetNpar()) {
563 
564  f1.ReleaseParameter(index);
565 
566  return true;
567 
568  } else {
569 
570  return false;
571  }
572  }
573 
574 
575  /**
576  * Set fit parameter limits.
577  *
578  * \param f1 fit function
579  * \param index parameter index
580  * \param xmin lower limit
581  * \param xmax upper limit
582  */
583  inline bool setParLimits(TF1& f1, const Int_t index, Double_t xmin, Double_t xmax)
584  {
585  using namespace std;
586 
587  if (index >= 0 && index < f1.GetNpar()) {
588 
589  if (xmin == 0.0) { xmin = -numeric_limits<Double_t>::min(); }
590  if (xmax == 0.0) { xmax = +numeric_limits<Double_t>::min(); }
591 
592  f1.SetParLimits(index, xmin, xmax);
593 
594  return true;
595 
596  } else {
597 
598  return false;
599  }
600  }
601 
602 
603  /**
604  * Check if fit parameter is fixed.
605  *
606  * \param f1 fit function
607  * \param index parameter index
608  */
609  inline bool isParameterFixed(const TF1& f1, const Int_t index)
610  {
611  if (index >= 0 && index < f1.GetNpar()) {
612 
613  Double_t xmin;
614  Double_t xmax;
615 
616  f1.GetParLimits(index, xmin, xmax);
617 
618  return (xmin != 0.0 && xmax != 0.0 && xmin >= xmax);
619 
620  } else {
621 
622  return false;
623  }
624  }
625 
626 
627  /**
628  * Get result of given textual formula.
629  *
630  * The formula may contain names of data members of the given object.
631  * For example:
632  * <pre>
633  * getResult("a / b.value", object, ..);
634  * </pre>
635  * In this, the corresponding data type should have (at least) data members <tt>a</tt> and <tt>b</tt>
636  * and <tt>b</tt> should have (at least) data member <tt>value</tt>.
637  *
638  * \param text text
639  * \param object object
640  * \param dictionary dictionary
641  * \return value
642  */
643  template<class T>
644  inline Double_t getResult(const TString& text,
645  const T& object,
646  const JRootDictionary_t& dictionary = JRootDictionary::getInstance())
647  {
648  using namespace std;
649  using namespace JPP;
650 
651  const TRegexp DOUBLE("[0-9.][eE][+-0-9]");
652 
653  string buffer = (const char*) text;
654 
655  for (string::size_type pos = 0; pos != buffer.size(); ) {
656 
657  if (isalpha(buffer[pos]) || buffer[pos] == '_') {
658 
659  string::size_type len = 1;
660 
661  while (pos + len != buffer.size() && (isalnum(buffer[pos+len]) || buffer[pos+len] == '_' || buffer[pos+len] == '.')) {
662  ++len;
663  }
664 
665  if (len != 1 || pos == 0 || TString(buffer.substr(pos-1).c_str()).Index(DOUBLE) != 0) {
666 
667  ostringstream os;
668 
669  os.setf(ios::fixed);
670  os.precision(10);
671 
672  JRootPrinter::print(os, object, buffer.substr(pos,len), dictionary);
673 
674  buffer.replace(pos, len, os.str());
675 
676  pos += os.str().size();
677 
678  continue;
679  }
680  }
681 
682  pos += 1;
683  }
684 
685  return TFormula("/tmp", buffer.c_str()).Eval(0.0);
686  }
687 
688 
689  /**
690  * Helper method to convert a 1D histogram to a graph.
691  *
692  * The graph consists of:
693  * - bin centers as x values;
694  * - bin contents as y values.
695  *
696  * Note that underflow and overflow bins are ignored.
697  *
698  * \param h1 1D histogram
699  * \return pointer to newly create graph
700  */
701  inline TGraph* histogramToGraph(const TH1& h1)
702  {
703  const int N = h1.GetNbinsX();
704 
705  std::vector<double> x(N), y(N);
706 
707  for (int i = 0; i < N; i++) {
708  x[i] = h1.GetBinCenter (i + 1);
709  y[i] = h1.GetBinContent(i + 1);
710  }
711 
712  return new TGraph(N, &x[0], &y[0]);
713  }
714 
715 
716  /**
717  * Helper method for ROOT histogram projections.
718  *
719  * \param h2 2D histogram
720  * \param xmin lower limit
721  * \param xmax upper limit
722  * \param projection projection (x|X|y|Y)
723  * \return pointer to newly created 1D histogram
724  */
725  inline TH1* projectHistogram(const TH2& h2, const Double_t xmin, const Double_t xmax, const char projection)
726  {
727  switch (projection) {
728 
729  case 'x':
730  case 'X':
731  return h2.ProjectionX("_px", h2.GetYaxis()->FindBin(xmin), h2.GetYaxis()->FindBin(xmax));
732 
733  case 'y':
734  case 'Y':
735  return h2.ProjectionY("_py", h2.GetXaxis()->FindBin(xmin), h2.GetXaxis()->FindBin(xmax));
736 
737  default:
738  return NULL;
739  }
740  }
741 }
742 
743 
744 /**
745  * Read regular expression from input stream
746  *
747  * \param in input stream
748  * \param object regular expression
749  * \return output stream
750  */
751 inline std::istream& operator>>(std::istream& in, TRegexp& object)
752 {
753  std::string buffer;
754 
755  if (in >> buffer) {
756  object = TRegexp(buffer.c_str());
757  }
758 
759  return in;
760 }
761 
762 
763 /**
764  * Write regular expression to output stream
765  *
766  * \param out output stream
767  * \param object regular expression
768  * \return output stream
769  */
770 inline std::ostream& operator<<(std::ostream& out, const TRegexp& object)
771 {
772  return out;
773 }
774 
775 #endif
776 
const double xmax
Definition: JQuadrature.cc:24
JFitParameter_t(const Int_t index, const Double_t value)
Constructor.
void AddPoint(TGraph *g1, const Double_t x, const Double_t y)
Add point to TGraph.
then usage $script[< detector identifier >< run range >]< QA/QCfile > nExample script to produce data quality plots nWhen a detector identifier and run range are data are downloaded from the database nand subsequently stored in the given QA QC file
Definition: JDataQuality.sh:19
char text[TEXT_SIZE]
Definition: elog.cc:72
bool releaseParameter(TF1 &f1, const Int_t index)
Release fit parameter.
Print objects in ASCII format using ROOT dictionary.
JFitParameter_t()
Default constructor.
then echo Enter input within $TIMEOUT_S seconds echo n User name
Definition: JCookie.sh:42
int getStreamerVersion(TFile *file, const char *const name)
Get ROOT streamer version of class with given name.
Definition: JRoot.hh:19
const TStreamerInfo * getStreamerInfo(TFile *file, const char *const name)
Get ROOT streamer information of class with given name.
Auxiliary class for a type holder.
Definition: JType.hh:19
bool resetObject(JManager< JKey_t, JValue_t > *object, const bool reset=false)
Reset JManager object.
Definition: JManager.hh:366
bool isParameterFixed(const TF1 &f1, const Int_t index)
Check if fit parameter is fixed.
boost::property_tree::ptree parse(std::string str)
Definition: configure.hh:24
const JPolynome f1(1.0, 2.0, 3.0)
Function.
const int n
Definition: JPolint.hh:786
Auxiliary data structure for a parameter index and its value.
TFile * getFile() const
Get file.
Definition: JRootFile.hh:66
I/O formatting auxiliaries.
Double_t getResult(const TString &text, TObject *object=NULL)
Get result of given textual formula.
bool setParameter(TF1 &f1, const JFitParameter_t &parameter)
Set fit parameter.
do set_variable OUTPUT_DIRECTORY $WORKDIR T
static void print(JRootWriter &writer, const JRootWritableClass &cls)
Write class contents to output.
Definition: JRootPrinter.hh:62
then awk string
ROOT input file.
Definition: JRootFile.hh:95
bool setParLimits(TF1 &f1, const Int_t index, Double_t xmin, Double_t xmax)
Set fit parameter limits.
TGraph * histogramToGraph(const TH1 &h1)
Helper method to convert a 1D histogram to a graph.
Type definition of ROOT based dictionary for ASCII I/O.
void reset(T &value)
Reset value.
const double xmin
Definition: JQuadrature.cc:23
std::istream & operator>>(std::istream &in, JAANET::JHead &header)
Read header from input.
Definition: JHead.hh:1829
then usage $script< input file >[option[primary[working directory]]] nWhere option can be N
Definition: JMuonPostfit.sh:40
static data_type & getInstance()
Get unique instance of template class.
Definition: JSingleton.hh:27
then fatal Wrong number of arguments fi set_variable DETECTOR $argv[1] set_variable INPUT_FILE $argv[2] eval JPrintDetector a $DETECTOR O IDENTIFIER eval JPrintDetector a $DETECTOR O SUMMARY JAcoustics sh $DETECTOR_ID source JAcousticsToolkit sh CHECK_EXIT_CODE typeset A EMITTERS get_tripods $WORKDIR tripod txt EMITTERS get_transmitters $WORKDIR transmitter txt EMITTERS for EMITTER in
Definition: JCanberra.sh:48
bool fixParameter(TF1 &f1, const JFitParameter_t &parameter)
Fix fit parameter.
std::ostream & operator<<(std::ostream &stream, const CLBCommonHeader &header)
const char * getName()
Get ROOT name of given data type.
Definition: JRootToolkit.hh:60
TH1 * projectHistogram(const TH2 &h2, const Double_t xmin, const Double_t xmax, const char projection)
Helper method for ROOT histogram projections.
Double_t g1(const Double_t x)
Function.
Definition: JQuantiles.cc:25