Jpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
JHalibut.cc File Reference

Example program to determine N-fold coincidence rates. More...

#include <string>
#include <iostream>
#include <fstream>
#include <iomanip>
#include <vector>
#include <map>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TF1.h"
#include "JDAQ/JDAQTimeslice.hh"
#include "JDAQ/JDAQClock.hh"
#include "JDAQ/JDAQEvaluator.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JModuleRouter.hh"
#include "JTrigger/JHitR0.hh"
#include "JTrigger/JBuildL0.hh"
#include "JTrigger/JSuperFrame1D.hh"
#include "JTrigger/JSuperFrame2D.hh"
#include "JGizmo/JManager.hh"
#include "JCalibrate/JCalibrateK40.hh"
#include "JROOT/JROOTClassSelector.hh"
#include "JTools/JQuantile.hh"
#include "JTools/JRange.hh"
#include "JSupport/JMultipleFileScanner.hh"
#include "JSupport/JTreeScanner.hh"
#include "JSupport/JAutoTreeScanner.hh"
#include "JSupport/JSupport.hh"
#include "JLang/JObjectMultiplexer.hh"
#include "JMath/JMathToolkit.hh"
#include "Jeep/JPrint.hh"
#include "Jeep/JParser.hh"
#include "Jeep/JMessage.hh"

Go to the source code of this file.

Functions

int main (int argc, char **argv)
 

Detailed Description

Example program to determine N-fold coincidence rates.

Author
mdejong

Definition in file JHalibut.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 51 of file JHalibut.cc.

52 {
53  using namespace std;
54  using namespace JPP;
55  using namespace KM3NETDAQ;
56 
57  typedef JRange<int> JRange_t;
58 
59  JMultipleFileScanner<> inputFile;
60  JLimit_t& numberOfEvents = inputFile.getLimit();
61  string outputFile;
62  string detectorFile;
63  double TMax_ns;
64  JROOTClassSelector selector;
65  JRange_t M;
66  string summaryFile;
67  string option;
68  int debug;
69 
70  try {
71 
72  JParser<> zap("Example program to determine N-fold coincidence rates.");
73 
74  zap['f'] = make_field(inputFile);
75  zap['o'] = make_field(outputFile) = "halibut.root";
76  zap['n'] = make_field(numberOfEvents) = JLimit::max();
77  zap['a'] = make_field(detectorFile);
78  zap['T'] = make_field(TMax_ns) = 20.0;
79  zap['C'] = make_field(selector) = getROOTClassSelection<JDAQTimesliceTypes_t>();
80  zap['M'] = make_field(M) = JRange_t(2, 31);
81  zap['s'] = make_field(summaryFile) = "halibut.txt";
82  zap['O'] = make_field(option) = "";
83  zap['d'] = make_field(debug) = 1;
84 
85  zap(argc, argv);
86  }
87  catch(const exception &error) {
88  FATAL(error.what() << endl);
89  }
90 
91  cout.tie(&cerr);
92 
93 
94  JDetector detector;
95 
96  try {
97  load(detectorFile, detector);
98  }
99  catch(const JException& error) {
100  FATAL(error);
101  }
102 
103  const JModuleRouter router(detector);
104 
105  const JMatch_t match(TMax_ns); // time window self-coincidences [ns]
106 
107  map<int, double> TMax_s; // histogram upper limit
108 
109  TMax_s[2] = 5.0e-2;
110  TMax_s[3] = 0.2e+0;
111  TMax_s[4] = 1.0e+0;
112  TMax_s[5] = 1.0e+1;
113  TMax_s[6] = 1.0e+2;
114 
115  typedef JSuperFrame2D<JHitR0> JSuperFrame2D_t;
116  typedef JSuperFrame1D<JHitR0> JSuperFrame1D_t;
117 
118  typedef JManager<int, TH1D> JManager_t;
119 
120  JManager_t H1(new TH1D("H1[%]", NULL, 100, -TMax_ns, +TMax_ns));
121  JManager_t T1(new TH1D("T1[%]", NULL, 100, 0.0, TMax_s[M.getLowerLimit()]));
122 
123  H1->Sumw2();
124  T1->Sumw2();
125 
127 
128  JAutoTreeScanner<JDAQTimeslice, JDAQEvaluator> zmap = JType<JDAQTimesliceTypes_t>();
129 
130  JTreeScannerInterface<JDAQTimeslice>* ps = zmap[selector];
131 
132  counter_type counter = 0;
133 
134  // process file-by-file to speed up JTreeScanner::configure();
135 
136  for (JMultipleFileScanner<>::const_iterator i = inputFile.begin(); i != inputFile.end(); ++i) {
137 
138  STATUS("Processing " << *i << endl);
139 
140  ps->configure(*i);
141 
142  vector<double> t0(detector.size(), 0.0);
143 
144  for ( ; ps->hasNext() && counter != inputFile.getLimit(); ++counter) {
145 
146  STATUS("event: " << setw(10) << counter << '\r'); DEBUG(endl);
147 
148  const JDAQTimeslice* timeslice = ps->next();
149 
150  for (JDAQTimeslice::const_iterator frame = timeslice->begin(); frame != timeslice->end(); ++frame) {
151 
152  JSuperFrame2D_t& buffer = JSuperFrame2D_t::demultiplex(*frame, router.getModule(frame->getModuleID()));
153 
154  for (JSuperFrame2D_t::iterator i = buffer.begin(); i != buffer.end(); ++i) {
155  i->join(match);
156  }
157 
158  JSuperFrame1D_t& data = JSuperFrame1D_t::multiplex(buffer);
159 
160  data.pop_back(); // remove end marker
161 
162  if (data.size() > 1) {
163 
164  TH1D* h1 = H1[frame->getModuleID()];
165  TH1D* t1 = T1[frame->getModuleID()];
166 
167  for (vector<JHitR0>::const_iterator p = data.begin(); p != data.end(); ) {
168 
170 
171  while (++q != data.end() && q->getT() - p->getT() <= TMax_ns ) {}
172 
173  const int N = distance(p,q);
174 
175  if (M(N)) {
176 
177  const int i = router.getIndex(frame->getModuleID());
178  const double ts = getTimeOfRTS(frame->getFrameIndex()) + p->getT();
179 
180  t1->Fill((ts - t0[i]) * 1.0e-9); // [s]
181 
182  t0[i] = ts;
183 
184  const double W = factorial(N, 2);
185 
186  for (vector<JHitR0>::const_iterator __p = p; __p != q; ++__p) {
187  for (vector<JHitR0>::const_iterator __q = __p; ++__q != q; ) {
188  h1->Fill(JCombinatorics::getSign(__p->getPMT(),__q->getPMT()) * (__q->getT() - __p->getT()), 1.0/W);
189  }
190  }
191  }
192 
193  p = q;
194  }
195  }
196  }
197  }
198  STATUS(endl);
199  }
200 
201  if (counter != 0) {
202 
203  const double V = (H1->GetXaxis()->GetXmax() - H1->GetXaxis()->GetXmin()) / (double) H1->GetXaxis()->GetNbins(); // [ns]
204  const double W = counter * getFrameTime() * 1.0e-9; // [s]
205 
206  for (JManager_t::iterator i = H1.begin(); i != H1.end(); ++i) {
207  i->second->Scale(1.0/(V*W));
208  }
209 
210  for (JManager_t::iterator i = T1.begin(); i != T1.end(); ++i) {
211  i->second->Scale(1.0/i->second->GetMaximum());
212  }
213  }
214 
215  if (summaryFile != "") {
216 
217  const double V = (H1->GetXaxis()->GetXmax() - H1->GetXaxis()->GetXmin()) / (double) H1->GetXaxis()->GetNbins(); // [ns]
218 
219  const int number_of_strings = getNumberOfStrings(detector);
220  const int number_of_floors = getNumberOfFloors (detector);
221  const int PRECISION = (M.getLowerLimit() > 2 ? 4 : 3);
222 
223  ofstream out(summaryFile.c_str());
224 
225  out << "Multiplicity " << M << endl;
226  out << "-------------------------------------------------------" << endl;
227  out << " location | Gauss | S - B | Total | slope " << endl;
228  out << " | [Hz] | [Hz] | [Hz] | [Hz] " << endl;
229  out << "-------------------------------------------------------" << endl;
230 
231  JQuantile Q[4];
232 
233  for (int string = 1; string <= number_of_strings; ++string) {
234  for (int floor = number_of_floors; floor >= 1; --floor) {
235 
236  const int id = detector.getModule(JModuleLocation(string,floor)).getID();
237 
238  out << " " << setw(3) << string << ' ' << setw(2) << floor << " ";
239 
240  TH1D* h1 = (H1.find(id) != H1.end() ? H1[id] : NULL);
241  TH1D* t1 = (T1.find(id) != T1.end() ? T1[id] : NULL);
242 
243  if (h1 != NULL) {
244 
245  TF1 f1("f1", "[0]*exp(-0.5*(x-[1])*(x-[1])/([2]*[2]))/(TMath::Sqrt(2*TMath::Pi())*[2]) + [3]");
246 
247  f1.SetParameter(0, h1->GetMaximum());
248  f1.SetParameter(1, 0.0);
249  f1.SetParameter(2, h1->GetRMS() * 0.25);
250  f1.SetParameter(3, h1->GetMinimum());
251 
252  h1->Fit(&f1, option.c_str(), "same");
253 
254  out << " | " << FIXED(8,PRECISION) << f1.GetParameter(0);
255  out << " | " << FIXED(8,PRECISION) << (h1->GetSumOfWeights() - f1.GetParameter(3) * h1->GetNbinsX()) * V;
256  out << " | " << FIXED(8,PRECISION) << h1->GetSumOfWeights() * V;
257 
258  Q[0].put( f1.GetParameter(0));
259  Q[1].put((h1->GetSumOfWeights() - f1.GetParameter(3) * h1->GetNbinsX()) * V);
260  Q[2].put( h1->GetSumOfWeights() * V);
261  }
262 
263  if (t1 != NULL) {
264 
265  TF1 f1("f1", "[0]*exp(-[1]*x)");
266 
267  f1.SetParameter(0, t1->GetMaximum());
268  f1.SetParameter(1, 1.0 / t1->GetRMS());
269 
270  t1->Fit(&f1, option.c_str(), "same");
271 
272  out << " | " << FIXED(8,PRECISION) << f1.GetParameter(1);
273 
274  Q[3].put(f1.GetParameter(1));
275  }
276 
277  out << endl;
278  }
279 
280  out << endl;
281  }
282 
283  if (Q[0].getCount() != 0) {
284 
285  out << "-------------------------------------------------------" << endl;
286  out << setw(10) << left << " average";
287 
288  for (int i = 0; i != sizeof(Q)/sizeof(Q[0]); ++i) {
289  out << " | " << FIXED(8,PRECISION) << Q[i].getMean();
290  }
291 
292  out << endl;
293  }
294 
295  out.close();
296  }
297 
298  if (outputFile != "") {
299 
300  TFile out(outputFile.c_str(), "RECREATE");
301 
302  H1.Write(out);
303  T1.Write(out);
304 
305  out.Close();
306  }
307 }
Utility class to parse command line options.
Definition: JParser.hh:1410
int getNumberOfStrings(const JDetector &detector)
Get number of strings.
double getCount(TH1D *hptr, int muon_threshold)
Auxiliary class to manage set of histograms.
#define STATUS(A)
Definition: JMessage.hh:61
Long64_t counter_type
Type definition for counter.
Definition: JCounter.hh:24
JSuperFrame2D< hit_type > JSuperFrame2D_t
Definition: JDataFilter.cc:93
Auxiliary class for a type holder.
Definition: JType.hh:19
Auxiliary data structure for floating point format specification.
Definition: JPrint.hh:461
double getTimeOfRTS(const JDAQChronometer &chronometer)
Get time of last RTS in ns since start of run for a given chronometer.
string outputFile
long long int factorial(const long long int n)
Determine factorial.
Definition: JMathToolkit.hh:40
JLimit JLimit_t
Type definition of limit.
Definition: JLimit.hh:214
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1836
double getFrameTime()
Get frame time duration.
Definition: JDAQClock.hh:162
int getNumberOfFloors(const JDetector &detector)
Get number of floors.
Data time slice.
void load(const JString &file_name, JDetector &detector)
Load detector from input file.
int debug
debug level
Definition: JSirene.cc:59
#define FATAL(A)
Definition: JMessage.hh:65
std::vector< frame_type >::iterator iterator
2-dimensional frame with time calibrated data from one optical module.
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:72
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:60
JSuperFrame1D< hit_type > JSuperFrame1D_t
Definition: JDataFilter.cc:92