Jpp
JHistogramToolkit.hh
Go to the documentation of this file.
1 #ifndef __JHISTOGRAM_TOOLKIT__
2 #define __JHISTOGRAM_TOOLKIT__
3 
4 #include <string>
5 
6 #include "TROOT.h"
7 #include "TFile.h"
8 #include "TH1D.h"
9 
10 #include "JFit/JFitParameters.hh"
11 
12 
13 /**
14  * \author mdejong
15  */
16 
17 using namespace JFIT;
18 
19 static const int JQUALITY = -1; //!< fit quality identifier
20 
21 
22 /**
23  * Map key.
24  */
25 struct JKey {
26  /**
27  * Constructor.
28  *
29  * \param name parameter name
30  * \param index parameter value (= index in user data)
31  */
32  JKey(const char* name,
33  const int index)
34  {
35  this->name = name;
36  this->index = index;
37  }
38 
39 
40  /**
41  * Less than operator.
42  *
43  * \param first first key
44  * \param second second key
45  * \return true if index of first key less than that of second; else false
46  */
47  friend inline bool operator<(const JKey& first, const JKey& second)
48  {
49  return first.index < second.index;
50  }
51 
52  std::string name;
53  int index;
54 };
55 
56 
57 /**
58  * Make key.
59  *
60  * \param PARAMETER parameter
61  * \return key
62  */
63 #define make_key(PARAMETER) JKey(#PARAMETER, PARAMETER)
64 
65 
66 /**
67  * Map value.
68  */
69 struct JValue {
70  /**
71  * Default constructor.
72  */
73  JValue() :
74  hA(NULL),
75  hB(NULL),
76  hC(NULL),
77  hD(NULL),
78  logx(false)
79  {}
80 
81 
82  /**
83  * Constructor.
84  *
85  * \param hA pointer to histogram
86  * \param hB pointer to histogram
87  * \param hC pointer to histogram
88  * \param hD pointer to histogram
89  * \param logx log10(x) filling mode
90  */
91  JValue(TH1D* hA,
92  TH1D* hB,
93  TH1D* hC,
94  TH1D* hD,
95  bool logx = false)
96  {
97  this->hA = hA;
98  this->hB = hB;
99  this->hC = hC;
100  this->hD = hD;
101  this->logx = logx;
102  }
103 
104 
105  /**
106  * Fill histograms.
107  *
108  * \param xA first abscissa value
109  * \param xB second abscissa value
110  * \param option option
111  */
112  void Fill(const Double_t xA, const Double_t xB, const bool option)
113  {
114  hA->Fill(logx ? log10(xA) : xA);
115  hB->Fill(logx ? log10(xB) : xB);
116  if (option)
117  hC->Fill(logx ? log10(xB/xA) : xB - xA);
118  else
119  hD->Fill(logx ? log10(xB/xA) : xB - xA);
120  }
121 
122 
123  /**
124  * Scale histograms.
125  */
126  void Scale()
127  {
128  if (hA->GetEntries() != 0) { hA->Scale(1.0/hA->GetEntries()); }
129  if (hB->GetEntries() != 0) { hB->Scale(1.0/hB->GetEntries()); }
130  if (hC->GetEntries() != 0) { hC->Scale(1.0/hC->GetEntries()); }
131  if (hD->GetEntries() != 0) { hD->Scale(1.0/hD->GetEntries()); }
132  }
133 
134 
135  /**
136  * Write histograms to file.
137  *
138  * \param out ROOT file
139  */
140  void Write(TFile& out)
141  {
142  out.WriteTObject(hA);
143  out.WriteTObject(hB);
144  out.WriteTObject(hC);
145  out.WriteTObject(hD);
146  }
147 
148 protected:
149  TH1D* hA;
150  TH1D* hB;
151  TH1D* hC;
152  TH1D* hD;
153  bool logx;
154 };
155 
156 
157 /**
158  * Auxiliary class to manage set of histograms.
159  */
160 struct JManager :
161  public std::map<JKey, JValue>
162 {
163  /**
164  * Insert set of histograms.
165  *
166  * \param key key
167  * \param nx number of bins
168  * \param xmin minimal abscissa value
169  * \param xmax maximal abscissa value
170  * \param logx log10(x)
171  */
172  void insert(const JKey& key,
173  const Int_t nx,
174  const Double_t xmin,
175  const Double_t xmax,
176  const double logx = false)
177  {
178  using namespace std;
179 
180  TH1D* hA = new TH1D(("[A " + key.name + "]").c_str(), NULL, nx, xmin, xmax);
181  TH1D* hB = new TH1D(("[B " + key.name + "]").c_str(), NULL, nx, xmin, xmax);
182  TH1D* hC = new TH1D(("[C " + key.name + "]").c_str(), NULL, nx, -xmax, xmax);
183  TH1D* hD = new TH1D(("[D " + key.name + "]").c_str(), NULL, nx, -xmax, xmax);
184 
185  std::map<JKey, JValue>::insert(make_pair(key, JValue(hA,hB,hC,hD,logx)));
186  }
187 
188 
189  /**
190  * Fill histograms.
191  *
192  * \param fA first fit result
193  * \param fB second fit result
194  * \param option option
195  */
196  void Fill(const JFit& fA, const JFit& fB, const bool option)
197  {
198  for (iterator i = begin(); i != end(); ++i) {
199 
200  const int index = i->first.index;
201 
202  if (index == JQUALITY) {
203  i->second.Fill(fA.getQ(), fB.getQ(), option);
204  } else if (fA.hasW(index) && fB.hasW(index)) {
205  i->second.Fill(fA.getW(index), fB.getW(index), option);
206  }
207  };
208  }
209 
210 
211  /**
212  * Scale histograms.
213  */
214  void Scale()
215  {
216  for (iterator i = this->begin(); i != this->end(); ++i) {
217  i->second.Scale();
218  }
219  }
220 
221 
222  /**
223  * Write histograms to file.
224  *
225  * \param out ROOT file
226  */
227  void Write(TFile& out)
228  {
229  for (iterator i = this->begin(); i != this->end(); ++i) {
230  i->second.Write(out);
231  }
232  }
233 
234 
235  /**
236  * Write histograms to file.
237  *
238  * \param file_name file name
239  */
240  void Write(const char* file_name)
241  {
242  TFile out(file_name, "RECREATE");
243 
244  this->Write(out);
245 
246  out.Write();
247  out.Close();
248  }
249 };
250 
251 #endif
JQUALITY
static const int JQUALITY
fit quality identifier
Definition: JHistogramToolkit.hh:19
JManager::Write
void Write(TFile &out)
Write histograms to file.
Definition: JHistogramToolkit.hh:227
std::iterator
Definition: JSTDTypes.hh:18
JKey::JKey
JKey(const char *name, const int index)
Constructor.
Definition: JHistogramToolkit.hh:32
JFIT
Auxiliary classes and methods for linear and iterative data regression.
Definition: JEnergy.hh:15
JManager::Scale
void Scale()
Scale histograms.
Definition: JHistogramToolkit.hh:214
JFIT::JFit::hasW
bool hasW(const int i) const
Check availability of value.
Definition: JEvt.hh:222
JValue::logx
bool logx
Definition: JHistogramToolkit.hh:153
JValue::hB
TH1D * hB
Definition: JHistogramToolkit.hh:150
JValue::Fill
void Fill(const Double_t xA, const Double_t xB, const bool option)
Fill histograms.
Definition: JHistogramToolkit.hh:112
JGIZMO::JManager
Auxiliary class to manage set of compatible ROOT objects (e.g.
Definition: JManager.hh:40
JManager::Write
void Write(const char *file_name)
Write histograms to file.
Definition: JHistogramToolkit.hh:240
JValue::Write
void Write(TFile &out)
Write histograms to file.
Definition: JHistogramToolkit.hh:140
JFIT::JFit::getW
const std::vector< double > & getW() const
Get values.
Definition: JEvt.hh:188
JFIT::JFit::getQ
double getQ() const
Get quality.
Definition: JEvt.hh:151
JValue::hC
TH1D * hC
Definition: JHistogramToolkit.hh:151
JValue::hA
TH1D * hA
Definition: JHistogramToolkit.hh:149
JKey::operator<
friend bool operator<(const JKey &first, const JKey &second)
Less than operator.
Definition: JHistogramToolkit.hh:47
JValue::hD
TH1D * hD
Definition: JHistogramToolkit.hh:152
JKey::name
std::string name
Definition: JHistogramToolkit.hh:52
std::map
Definition: JSTDTypes.hh:16
JGIZMO::JManager::Write
void Write(TDirectory &out, const bool wm=false)
Write objects to file.
Definition: JManager.hh:258
JKey::index
int index
Definition: JHistogramToolkit.hh:53
JValue::JValue
JValue()
Default constructor.
Definition: JHistogramToolkit.hh:73
JManager::Fill
void Fill(const JFit &fA, const JFit &fB, const bool option)
Fill histograms.
Definition: JHistogramToolkit.hh:196
JFitParameters.hh
std
Definition: jaanetDictionary.h:36
JValue::Scale
void Scale()
Scale histograms.
Definition: JHistogramToolkit.hh:126
JValue::JValue
JValue(TH1D *hA, TH1D *hB, TH1D *hC, TH1D *hD, bool logx=false)
Constructor.
Definition: JHistogramToolkit.hh:91
JManager::insert
void insert(const JKey &key, const Int_t nx, const Double_t xmin, const Double_t xmax, const double logx=false)
Insert set of histograms.
Definition: JHistogramToolkit.hh:172
JKey
Map key.
Definition: JHistogramToolkit.hh:25
JFIT::JFit
Data structure for track fit results.
Definition: JEvt.hh:31
JLANG::JValue
Wrapper class around template object.
Definition: JValue.hh:151