Jpp  15.0.3
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JMergeRunAnalyzer.cc
Go to the documentation of this file.
1 #include <string>
2 #include <iostream>
3 #include <iomanip>
4 #include <vector>
5 #include <map>
6 
7 #include "TROOT.h"
8 #include "TFile.h"
9 #include "TH1D.h"
10 #include "TH2D.h"
11 #include "TKey.h"
12 #include "TString.h"
13 #include "TMath.h"
14 #include "TTimeStamp.h"
15 
17 
18 #include "JROOT/JRootFileReader.hh"
19 #include "JROOT/JRootFileWriter.hh"
20 #include "JSupport/JMeta.hh"
21 
23 #include "JTrigger/JTriggerBits.hh"
24 
25 #include "JGizmo/JGizmoToolkit.hh"
26 
27 #include "Jeep/JParser.hh"
28 #include "Jeep/JMessage.hh"
29 
30 /**
31  * \author adomi
32  * Auxiliary program to merge JRunAnalyzer histograms.
33  */
34 
35 /*
36  * Fill TH2 from TH1 of each run.
37  */
38 inline void fill2DHistogram(int run, TH1D* h1, TH2D* h2, bool norm=false){
39 
40  double normalisation = 1;
41  if(norm) normalisation = h1->Integral();
42 
43  for(int i = 1; i <= h1->GetNbinsX(); ++i){
44 
45  h2->Fill(run, h1->GetXaxis()->GetBinCenter(i), h1->GetBinContent(i)/normalisation);
46  }
47 }
48 
49 int main(int argc, char **argv)
50 {
51  using namespace std;
52  using namespace JPP;
53  using namespace KM3NETDAQ;
54 
55  vector<string> inputFile;
56  string outputFile;
57  int debug;
58 
59  try {
60 
61  JParser<> zap("Auxiliary program to merge JRunAnalyzer histograms.");
62 
63  zap['f'] = make_field(inputFile, "input file (output from JRunAnalyzer).");
64  zap['o'] = make_field(outputFile, "output file.") = "merge-jra.root";
65  zap['d'] = make_field(debug, "debug.") = 1;
66 
67  zap(argc, argv);
68  }
69  catch(const exception &error) {
70  FATAL(error.what() << endl);
71  }
72 
73  int min_run = numeric_limits<int>::max(), max_run = 0;
74 
75  size_t min_date = numeric_limits<int>::max(), max_date = 0;
76 
77  for (vector<string>::const_iterator i = inputFile.begin(); i != inputFile.end(); ++i) {
78 
79  TFile in(i->c_str(), "read");
80 
81  for (JRootFileReader<JDAQHeader> in(i->c_str()); in.hasNext(); ) {
82 
83  const JDAQHeader* p = in.next();
84 
85  if(p->getTimesliceStart().getUTCseconds() < min_date) min_date = p->getTimesliceStart().getUTCseconds();
86  if(p->getTimesliceStart().getUTCseconds() > max_date) max_date = p->getTimesliceStart().getUTCseconds();
87  if(p->getRunNumber() < min_run) min_run = p->getRunNumber();
88  if(p->getRunNumber() > max_run) max_run = p->getRunNumber();
89 
90  }
91 
92  in.Close();
93  }
94 
95  time_t tmin = min_date, tmax = max_date;
96 
97  TTimeStamp date_min(tmin, 0);
98  TTimeStamp date_max(tmax, 0);
99 
100  cout << "Minimum RUN: " << min_run << ", "; date_min.Print("s");
101  cout << "Maximum RUN: " << max_run << ", "; date_max.Print("s");
102 
103  TFile out(outputFile.c_str(), "recreate");
104 
105  int xbins = (max_run - min_run) + 1;
106  string title = string("Data taken from: ") + date_min.AsString("s") + " - " + date_max.AsString("s");
107 
108  TH2D* h_pmt_rate_distribution = new TH2D("h_pmt_rate_distribution", string(title + "; RUN; PMT Rate Distribution [kHz]").c_str(), xbins, (double) min_run-0.5, (double) max_run+0.5, 20, 0, log10(100));
109  setLogarithmicY(h_pmt_rate_distribution);
110  h_pmt_rate_distribution -> SetMinimum(1);
111 
112  TH2D* h_Trigger_bit_hit = new TH2D("h_Trigger_bit_hit", string(title + "; RUN; Hits trigger bit").c_str(), xbins, (double) min_run-0.5, (double) max_run+0.5, NUMBER_OF_TRIGGER_BITS, -0.5, NUMBER_OF_TRIGGER_BITS - 0.5);
113 
114  TH2D* h_Triggered_hits = new TH2D("h_Triggered_hits", string(title + "; RUN; Number of triggered hits").c_str(), xbins, (double) min_run-0.5, (double) max_run+0.5, 50, 0, 4);
115  setLogarithmicY(h_Triggered_hits);
116 
117  TH2D* h_Triggered_hits_3dmuon = new TH2D("h_Triggered_hits_3dmuon", string(title + "; RUN; Number of triggered hits - JTRIGGER3DMUON").c_str(), xbins, (double) min_run-0.5, (double) max_run+0.5, 50, 0, 3);
118  setLogarithmicY(h_Triggered_hits_3dmuon);
119 
120  TH2D* h_Snapshot_hits = new TH2D("h_Snapshot_hits", string(title + "; RUN; Number of snapshot hits").c_str(), xbins, (double) min_run-0.5, (double) max_run+0.5, 50, 0, 4);
121  setLogarithmicY(h_Snapshot_hits);
122 
123  TH2D* h_Triggered_over_Snapshot_hits = new TH2D("h_Triggered_over_Snapshot_hits", string(title + "; RUN; Triggered/Snapshot hits").c_str(), xbins, (double) min_run-0.5, (double) max_run+0.5, 100, 0, 0.5);
124 
125  TH2D* h_Number_of_overlays = new TH2D("h_Number_of_overlays", string(title + "; RUN; Number of Overlays").c_str(), xbins, (double) min_run-0.5, (double) max_run+0.5, 1000, -0.5, 1000 - 0.5);
126 
127  TH2D* h_pmt_distribution_triggered_hits = new TH2D("h_pmt_distribution_triggered_hits", string(title + "; RUN;Normalised PMT Distrib triggered hits (upper PMTs=[0-11])").c_str(), xbins, (double) min_run-0.5, (double) max_run+0.5, NUMBER_OF_PMTS, -0.5 , NUMBER_OF_PMTS - 0.5);
128 
129  for (vector<string>::const_iterator i = inputFile.begin(); i != inputFile.end(); ++i) {
130 
131  DEBUG("Processing " << *i << endl) ;
132 
133  TFile in(i->c_str(), "read");
134 
135  int run = 0;
136 
137  for (JRootFileReader<JDAQHeader> in(i->c_str()); in.hasNext(); ) {
138 
139  const JDAQHeader* p = in.next();
140 
141  run = p->getRunNumber();
142 
143  }
144 
145  TH1D* h1d_prd = new TH1D();
146  TH1D* h1d_thb = new TH1D();
147  TH1D* h1d_th = new TH1D();
148  TH1D* h1d_th3D = new TH1D();
149  TH1D* h1d_sh = new TH1D();
150  TH1D* h1d_tosh = new TH1D();
151  TH1D* h1d_noo = new TH1D();
152  TH1D* h1d_pdth = new TH1D();
153 
154  in.GetDirectory("Detector")->GetObject("h_pmt_rate_distribution", h1d_prd);
155  fill2DHistogram(run, h1d_prd, h_pmt_rate_distribution);
156 
157  TDirectory* dir = in.GetDirectory("JDAQEvent");
158 
159  dir->GetObject("h_Trigger_bit_hit", h1d_thb);
160  fill2DHistogram(run, h1d_thb, h_Trigger_bit_hit);
161 
162  dir->GetObject("h_Triggered_hits", h1d_th);
163  fill2DHistogram(run, h1d_th, h_Triggered_hits);
164 
165  dir->GetObject("h_Triggered_hits_3dmuon", h1d_th3D);
166  fill2DHistogram(run, h1d_th3D, h_Triggered_hits_3dmuon);
167 
168  dir->GetObject("h_Snapshot_hits", h1d_sh);
169  fill2DHistogram(run, h1d_sh, h_Snapshot_hits);
170 
171  dir->GetObject("h_Triggered_over_Snapshot_hits", h1d_tosh);
172  fill2DHistogram(run, h1d_tosh, h_Triggered_over_Snapshot_hits);
173 
174  dir->GetObject("h_Number_of_overlays", h1d_noo);
175  fill2DHistogram(run, h1d_noo, h_Number_of_overlays);
176 
177  dir->GetObject("h_pmt_distribution_triggered_hits", h1d_pdth);
178  fill2DHistogram(run, h1d_pdth, h_pmt_distribution_triggered_hits, true);
179 
180  in.Close();
181 
182  }
183 
184  out.Write();
185  out.Close();
186 
187 }
Utility class to parse command line options.
Definition: JParser.hh:1500
int main(int argc, char *argv[])
Definition: Main.cc:15
JDAQUTCExtended getTimesliceStart() const
Get start of timeslice.
static const unsigned int NUMBER_OF_TRIGGER_BITS
Number of trigger bits.
then for HISTOGRAM in h0 h1
Definition: JMatrixNZ.sh:71
string outputFile
int getRunNumber() const
Get run number.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1961
ROOT I/O of application specific meta data.
int debug
debug level
Definition: JSirene.cc:63
JUINT32_t getUTCseconds() const
Get time.
General purpose messaging.
#define FATAL(A)
Definition: JMessage.hh:67
void fill2DHistogram(int run, TH1D *h1, TH2D *h2, bool norm=false)
Utility class to parse command line options.
then usage $script< input_file >< detector_file > fi set_variable OUTPUT_DIR set_variable SELECTOR JDAQTimesliceL1 set_variable DEBUG case set_variable DEBUG
void setLogarithmicY(TF2 *f2)
Make parameter y of function logarithmic (e.g. after filling with log10()).
Setting of trigger bits.
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition: JDAQ.hh:26
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 source JAcoustics sh $DETECTOR_ID CHECK_EXIT_CODE typeset A TRIPODS get_tripods $WORKDIR tripod txt TRIPODS for EMITTER in
Definition: JCanberra.sh:42
ROOT file reader.