Jpp  18.3.1
the software that should make you happy
 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/JDAQTimesliceIO.hh"
#include "km3net-dataformat/online/JDAQClock.hh"
#include "JDAQ/JDAQEvaluator.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JModuleRouter.hh"
#include "JTrigger/JHitR0.hh"
#include "JTrigger/JMatchL0.hh"
#include "JTrigger/JBuildL0.hh"
#include "JTrigger/JSuperFrame1D.hh"
#include "JTrigger/JSuperFrame2D.hh"
#include "JROOT/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 52 of file JHalibut.cc.

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