Jpp  19.1.0-rc.1
the software that should make you happy
Functions
JSquid.cc File Reference

Example program to histogram string and floor time difference. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <map>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "JTools/JQuantile.hh"
#include "JTools/JRange.hh"
#include "JROOT/JManager.hh"
#include "JROOT/JRootToolkit.hh"
#include "JDAQ/JDAQEventIO.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JDAQHitRouter.hh"
#include "JPhysics/JConstants.hh"
#include "JSupport/JMultipleFileScanner.hh"
#include "JSupport/JSupport.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 histogram string and floor time difference.

Author
mdejong

Definition in file JSquid.cc.

Function Documentation

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 36 of file JSquid.cc.

37 {
38  using namespace std;
39  using namespace JPP;
40  using namespace KM3NETDAQ;
41 
42  typedef JRange<double> JRange_t;
43 
45  JLimit_t& numberOfEvents = inputFile.getLimit();
46  string outputFile;
47  string detectorFile;
48  JRange_t T_ns;
49  double Wmin;
50  double Qmin;
51  int qaqc;
52  int debug;
53 
54  try {
55 
56  JParser<> zap("Example program to histogram string and floor time difference.");
57 
58  zap['f'] = make_field(inputFile);
59  zap['o'] = make_field(outputFile) = "squid.root";
60  zap['a'] = make_field(detectorFile);
61  zap['n'] = make_field(numberOfEvents) = JLimit::max();
62  zap['T'] = make_field(T_ns, "Time window for coincidences [ns]") = JRange_t(-50.0, +200.0);
63  zap['W'] = make_field(Wmin, "Minimal number of entries") = 10.0;
64  zap['q'] = make_field(Qmin, "Minimal fraction of coincidences") = 0.5;
65  zap['Q'] = make_field(qaqc) = 0;
66  zap['d'] = make_field(debug) = 2;
67 
68  zap(argc, argv);
69  }
70  catch(const exception& error) {
71  FATAL(error.what() << endl);
72  }
73 
74 
76 
77  try {
78  load(detectorFile, detector);
79  }
80  catch(const JException& error) {
81  FATAL(error);
82  }
83 
84  const JDAQHitRouter router(detector);
85 
86  JManager<JLocation, TH1D> H1(new TH1D("%", NULL, 200, -1.0e3, +1.0e3));
87 
88  while (inputFile.hasNext()) {
89 
90  STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl);
91 
92  JDAQEvent* event = inputFile.next();
93 
95 
96  for (JDAQEvent::const_iterator<JDAQTriggeredHit> hit = event->begin<JDAQTriggeredHit>(); hit != event->end<JDAQTriggeredHit>(); ++hit) {
97  buffer[router.getModule(*hit)].put(getTime(*hit, router.getPMT(*hit)) + router.getModule(*hit).getZ() * getInverseSpeedOfLight());
98  }
99 
100  if (buffer.size() >= 2u) {
101  for (map<JLocation, JQuantile>::const_iterator q = buffer.begin(), p = q++; q != buffer.end(); ++p, ++q) {
102  if (p->first.getString() == q->first.getString() && p->first.getFloor() + 1 == q->first.getFloor()) {
103  H1[p->first]->Fill(q->second.getMean() - p->second.getMean());
104  }
105  }
106  }
107  }
108  STATUS(endl);
109 
110 
111  TFile out(outputFile.c_str(), "recreate");
112 
113  out << H1;
114 
115  out.Write();
116  out.Close();
117 
118 
119  int nin = 0;
120  int nout = 0;
121 
122  for (const auto& i : H1) {
123 
124  Double_t W = 0.0;
125 
126  TH1D* h1 = i.second;
127 
128  if (h1->GetSumOfWeights() > Wmin) {
129 
130  for (Int_t ix = 1; ix <= h1->GetXaxis()->GetNbins(); ++ix) {
131 
132  const Double_t x = h1->GetBinCenter (ix);
133  const Double_t y = h1->GetBinContent(ix);
134 
135  if (T_ns(x)) {
136  W += y;
137  }
138  }
139 
140  DEBUG(i.first << ' ' << FIXED(6,0) << W << '/' << FIXED(6,0) << h1->GetSumOfWeights() << endl);
141 
142  if (W / h1->GetSumOfWeights() >= Qmin)
143  nin += 1;
144  else
145  nout += 1;
146  }
147  }
148 
149  NOTICE("Number of modules out/in micro-sync " << nout << '/' << nin << endl);
150 
151  QAQC(nin << ' ' << nout << endl);
152 
153  return 0;
154 }
string outputFile
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62
#define STATUS(A)
Definition: JMessage.hh:63
#define QAQC(A)
QA/QC output macro.
Definition: JMessage.hh:100
int qaqc
QA/QC file descriptor.
#define NOTICE(A)
Definition: JMessage.hh:64
#define FATAL(A)
Definition: JMessage.hh:67
int debug
debug level
Definition: JSirene.cc:69
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:2158
Simple wrapper around JModuleRouter class for direct addressing of PMT data in detector data structur...
Detector data structure.
Definition: JDetector.hh:96
General exception.
Definition: JException.hh:24
Utility class to parse command line options.
Definition: JParser.hh:1714
Auxiliary class to manage set of compatible ROOT objects (e.g. histograms) using unique keys.
Definition: JManager.hh:47
General purpose class for object reading from a list of file names.
virtual bool hasNext() override
Check availability of next element.
counter_type getCounter() const
Get counter.
virtual const pointer_type & next() override
Get next element.
Template const_iterator.
Definition: JDAQEvent.hh:68
double getTime(const Hit &hit)
Get true time of hit.
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
const double getInverseSpeedOfLight()
Get inverse speed of light.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
KM3NeT DAQ data structures and auxiliaries.
Definition: DataQueue.cc:39
Definition: JSTDTypes.hh:14
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:448
Type definition of range.
Definition: JHead.hh:43
Detector file.
Definition: JHead.hh:227
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:45