Jpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
JMergeCalibrateK40.cc File Reference
#include <string>
#include <iostream>
#include <iomanip>
#include <vector>
#include <map>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TH2D.h"
#include "TKey.h"
#include "TString.h"
#include "TMath.h"
#include "JROOT/JRootFileReader.hh"
#include "JSupport/JMeta.hh"
#include "JCalibrate/JCalibrateK40.hh"
#include "Jeep/JParser.hh"
#include "Jeep/JMessage.hh"

Go to the source code of this file.

Functions

int main (int argc, char **argv)
 

Function Documentation

int main ( int  argc,
char **  argv 
)
Author
mkarel Auxiliary program to merge K40 data.

Definition at line 28 of file JMergeCalibrateK40.cc.

29 {
30  using namespace std;
31  using namespace JPP;
32 
33  vector<string> inputFile;
34  string outputFile;
35  int debug;
36 
37  try {
38 
39  JParser<> zap("Auxiliary program to merge K40 data.");
40 
41  zap['f'] = make_field(inputFile, "input file (output from JCalibrateK40).");
42  zap['o'] = make_field(outputFile, "output file (input to JFitK40).") = "mergek40.root";
43  zap['d'] = make_field(debug, "debug.") = 1;
44 
45  zap(argc, argv);
46  }
47  catch(const exception &error) {
48  FATAL(error.what() << endl);
49  }
50 
51 
52  TH1* weights_hist = NULL;
53 
54  typedef map<TString, TH1*> map_type;
55 
56  map_type zmap;
57 
58 
59  for (vector<string>::const_iterator i = inputFile.begin(); i != inputFile.end(); ++i) {
60 
61  DEBUG("Processing " << *i << endl) ;
62 
63  TFile in(i->c_str(), "read");
64 
65  // check if file contains weights vector
66 
67  if (!in.GetListOfKeys()->Contains(weights_hist_t)) {
68  FATAL("Histogram does not contain weights histogram." << endl);
69  }
70 
71  TH1* weights_hist_i = (TH1*) in.Get(weights_hist_t);
72 
73  if (weights_hist == NULL) {
74 
75  weights_hist = (TH1*) weights_hist_i->Clone();
76 
77  } else {
78 
79  // check normalization
80 
81  TAxis* x_axis = weights_hist->GetXaxis();
82 
83  if (weights_hist_i->GetBinContent(x_axis->FindBin(WS_t)) !=
84  weights_hist ->GetBinContent(x_axis->FindBin(WS_t))) {
85  FATAL("Number of y-bins does not match for this file " << *i << endl);
86  }
87 
88  if (weights_hist_i->GetBinContent(x_axis->FindBin(WB_t)) !=
89  weights_hist ->GetBinContent(x_axis->FindBin(WB_t))) {
90  FATAL("TMax_ns is not the same for this file " << *i << endl);
91  }
92 
93  // add normalization
94 
95  weights_hist->SetBinContent(x_axis->FindBin(W1_overall_t),
96  weights_hist ->GetBinContent(x_axis->FindBin(W1_overall_t)) +
97  weights_hist_i->GetBinContent(x_axis->FindBin(W1_overall_t)));
98  }
99 
100 
101  TIter iter(in.GetListOfKeys());
102 
103  for (TKey* key; (key = (TKey*) iter.Next()) != NULL; ) {
104 
105  if (TString(key->GetName()) != TString(weights_hist->GetName())) {
106 
107  TH1* h1 = dynamic_cast<TH1*>(key->ReadObj());
108 
109  if (h1 != NULL) {
110 
111  map_type::iterator p = zmap.find(h1->GetName());
112 
113  if (p == zmap.end()) {
114 
115  DEBUG("Clone " << h1->GetName() << endl);
116 
117  p = zmap.insert(make_pair(h1->GetName(), (TH1*) h1->Clone())).first;
118 
119  } else {
120 
121  DEBUG("Add " << h1->GetName() << endl);
122 
123  p->second->Add(h1);
124  }
125  }
126  }
127  }
128 
129  weights_hist->SetDirectory(0);
130 
131  for (map_type::iterator i = zmap.begin(); i != zmap.end(); ++i) {
132  i->second->SetDirectory(0);
133  }
134 
135  in.Close();
136  }
137 
138 
139  if (weights_hist == NULL) {
140  FATAL("Missing histogram " << weights_hist_t << endl);
141  }
142 
143 
144  // write file
145 
146  TFile out(outputFile.c_str(), "recreate");
147 
148  const Double_t WS = weights_hist->GetBinContent(weights_hist->GetXaxis()->FindBin(WS_t)) ; // [ns]
149  const Double_t WB = weights_hist->GetBinContent(weights_hist->GetXaxis()->FindBin(WB_t)) ; // [ns]
150 
151  for (map_type::iterator i = zmap.begin(); i != zmap.end(); ++i) {
152 
153  if (i->first.EndsWith(_2S)) {
154 
155  TH2D* h2s = (TH2D*) i->second;
156 
157  if (h2s->GetEntries() != 0) {
158 
159  TH1D* h1b = (TH1D*) zmap.find(TString(i->first).ReplaceAll(_2S, _1B))->second;
160  TH1D* h1L = (TH1D*) zmap.find(TString(i->first).ReplaceAll(_2S, _1L))->second;
161 
162  //----------------------------------------------------------
163  // normalisation and background estimation
164  //----------------------------------------------------------
165 
166  for (int ix = 1; ix <= h2s->GetXaxis()->GetNbins(); ++ix) {
167 
168  // normalised background estimation
169  // if live time is zero then set bin count to zero and error to some large value
170 
171  if (h1L->GetBinContent(ix) == 0.0) {
172 
173  h1b->SetBinContent(ix, 0.0);
174  h1b->SetBinError (ix, 0.0);
175 
176  for (int iy = 1; iy <= h2s->GetYaxis()->GetNbins(); ++iy) {
177  h2s->SetBinContent(ix, iy, 0.0);
178  h2s->SetBinError (ix, iy, 999.999);
179  }
180 
181  } else {
182 
183  const double Y = h1b->GetBinContent(ix) / h1L->GetBinContent(ix) / WB;
184  const double W = h1b->GetBinError (ix) / h1L->GetBinContent(ix) / WB;
185 
186  h1b->SetBinContent(ix, Y);
187  h1b->SetBinError (ix, W);
188 
189  for (int iy = 1; iy <= h2s->GetYaxis()->GetNbins(); ++iy) {
190 
191  double y = h2s->GetBinContent(ix, iy) / h1L->GetBinContent(ix) / WS;
192  double w = h2s->GetBinError (ix, iy) / h1L->GetBinContent(ix) / WS;
193 
194  if (w == 0.0) {
195 
196  // for chi2 minimisation; set accordingly error of empty bins
197 
198  w = 1.0 / h1L->GetBinContent(ix) / WS;
199  }
200 
201  h2s->SetBinContent(ix, iy, y - Y);
202  h2s->SetBinError (ix, iy, TMath::Sqrt(w*w + W*W));
203  }
204  }
205  }
206 
207  // reset under and overflows.
208 
209  for (int ix = 0; ix <= h2s->GetXaxis()->GetNbins() + 1; ++ix) {
210  h2s->SetBinContent(ix, 0, 0.0); h2s->SetBinContent(ix, h2s->GetYaxis()->GetNbins() + 1, 0.0);
211  h2s->SetBinError (ix, 0, 0.0); h2s->SetBinError (ix, h2s->GetYaxis()->GetNbins() + 1, 0.0);
212  }
213 
214  for (int iy = 0; iy <= h2s->GetYaxis()->GetNbins() + 1; ++iy) {
215  h2s->SetBinContent(0, iy, 0.0); h2s->SetBinContent(h2s->GetXaxis()->GetNbins() + 1, iy, 0.0);
216  h2s->SetBinError (0, iy, 0.0); h2s->SetBinError (h2s->GetXaxis()->GetNbins() + 1, iy, 0.0);
217  }
218 
219  NOTICE("Biased coincidence rate [Hz] " << h2s->GetName() << ' ' << h2s->GetSumOfWeights() * WS << endl);
220 
221  h2s->SetName(TString(i->first).ReplaceAll(_2S, _2R));
222 
223  h2s->Write();
224  }
225  }
226  }
227 
228  putObject(&out, JMeta(argc, argv));
229 
230  for (vector<string>::const_iterator i = inputFile.begin(); i != inputFile.end(); ++i) {
231  JMeta::copy(i->c_str(), out);
232  }
233 
234  out.Write();
235  out.Close();
236 }
Utility class to parse command line options.
Definition: JParser.hh:1410
static const char *const WB_t
Named bin for value of TMax_ns in JCalibrateK40.
string outputFile
static const char *const WS_t
Named bin for time residual bin width.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1836
static const char *const _1B
Name extension for 1D background.
#define NOTICE(A)
Definition: JMessage.hh:62
static const char *const _2R
Name extension for 2D rate.
int debug
debug level
Definition: JSirene.cc:59
#define FATAL(A)
Definition: JMessage.hh:65
static const char *const weights_hist_t
Histogram naming.
static const char *const _1L
Name extension for 1D live time.
static const char *const W1_overall_t
Named bin for duration of the run.
void copy(const Head &from, JHead &to)
Copy header from from to to.
Definition: JHead.cc:40
static const char *const _2S
Name extension for 2D counts.
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:60
bool putObject(TDirectory *dir, const T &object)
Write object to ROOT directory.