Jpp  master_rocky
the software that should make you happy
JDiffToA.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <iomanip>
3 #include <vector>
4 #include <algorithm>
5 
6 #include "TROOT.h"
7 #include "TFile.h"
8 
9 #include "JAcoustics/JToA.hh"
10 #include "JAcoustics/JSupport.hh"
11 
13 
14 #include "Jeep/JParser.hh"
15 #include "Jeep/JMessage.hh"
16 
17 namespace {
18 
19  using JACOUSTICS::JToA;
20 
21  /**
22  * Compare measured times of arrival.
23  *
24  * \param first first time-of-arrival
25  * \param second second time-of-arrival
26  * \return true if first time-of-arrival less than second; else false
27  */
28  inline bool compare(const JToA& first, const JToA& second)
29  {
30  if (first.DETID == second.DETID) {
31 
32  if (first.RUN == second.RUN) {
33 
34  if (first.DOMID == second.DOMID) {
35 
36  if (first.WAVEFORMID == second.WAVEFORMID) {
37 
38  if (first.TOA_NS == second.TOA_NS) {
39 
40  if (first.SECONDS == second.SECONDS) {
41 
42  if (first.TICKS == second.TICKS) {
43 
44  if (first.QUALITYFACTOR == second.QUALITYFACTOR ) {
45 
46  return first.QUALITYNORMALISATION < second.QUALITYNORMALISATION;
47  }
48  else {
49 
50  return first.QUALITYFACTOR < second.QUALITYFACTOR;
51  }
52  }
53  else {
54 
55  return first.TICKS < second.TICKS;
56  }
57  } else {
58 
59  return first.SECONDS < second.SECONDS;
60  }
61 
62  } else {
63 
64  return first.TOA_NS < second.TOA_NS;
65  }
66 
67  } else {
68 
69  return first.WAVEFORMID < second.WAVEFORMID;
70  }
71 
72  } else {
73 
74  return first.DOMID < second.DOMID;
75  }
76 
77  } else {
78 
79  return first.RUN < second.RUN;
80  }
81 
82  } else {
83 
84  return first.DETID < second.DETID;
85  }
86  }
87 
88  /**
89  * Write time-of-arrival to outout stream.
90  *
91  * \param out output stream
92  * \param object time-of-arrival
93  * \return output stream
94  */
95  std::ostream& operator<<(std::ostream& out, const JToA& object)
96  {
97  using namespace std;
98  using namespace JPP;
99 
100  out << setw(8) << object.RUN << ' '
101  << setw(8) << object.DOMID << ' '
102  << setw(2) << object.WAVEFORMID << ' '
103  << FIXED(15,1) << (static_cast<double>(object.SECONDS) + static_cast<double>(object.TICKS)*16E-9) << ' '
104  << FIXED( 9,6) << object.TOA_NS << ' '
105  << FIXED( 9,3) << object.QUALITYFACTOR;
106 
107  return out;
108  }
109 }
110 
111 /**
112  * \file
113  *
114  * Program to compare toa data.
115  * \author mdejong
116  */
117 int main(int argc, char **argv)
118 {
119  using namespace std;
120  using namespace JPP;
121 
122  vector<string> inputFile;
123  JLimit numberOfEvents;
124  int debug;
125 
126  try {
127 
128  JParser<> zap("Program to compare toa data.");
129 
130  zap['f'] = make_field(inputFile, "two outputs of JToA");
131  zap['n'] = make_field(numberOfEvents) = JLimit::max();
132  zap['d'] = make_field(debug) = 2;
133 
134  zap(argc, argv);
135  }
136  catch(const exception &error) {
137  FATAL(error.what() << endl);
138  }
139 
140  if (inputFile.size() != 2u) {
141  FATAL("Wrong number of input files " << inputFile.size() << endl);
142  }
143 
144  const size_t width = max(inputFile[0].size(), inputFile[1].size());
145 
146  vector<JToA> buffer[2];
147 
148  for (int i = 0; i != 2; ++i) {
149 
150  for (JSingleFileScanner<JToA> in(inputFile[i], numberOfEvents); in.hasNext(); ) {
151  buffer[i].push_back(*in.next());
152  }
153 
154  sort(buffer[i].begin(), buffer[i].end(), compare);
155  }
156 
157  int count[] = { 0, 0 };
158 
160  p0 = buffer[0].begin(),
161  p1 = buffer[1].begin(); p0 != buffer[0].end() && p1 != buffer[1].end(); ) {
162 
163  for ( ; p0 != buffer[0].end() && p1 != buffer[1].end() && compare(*p0,*p1); ++p0, ++count[1]) {
164  DEBUG(">> " << setw(width) << left << inputFile[0] << right << ' ' << *p0 << endl);
165  }
166 
167  for ( ; p0 != buffer[0].end() && p1 != buffer[1].end() && compare(*p1,*p0); ++p1, ++count[1]) {
168  DEBUG("<< " << setw(width) << left << inputFile[1] << right << ' ' << *p1 << endl);
169  }
170 
171  if (p0 != buffer[0].end() && p1 != buffer[1].end() && !compare(*p0,*p1) && !compare(*p1,*p0)) {
172 
173  ++count[0];
174 
175  DEBUG(setw(width) << left << inputFile[0] << right << ' ' << *p0 << " \\" << endl);
176  DEBUG(setw(width) << left << inputFile[1] << right << ' ' << *p1 << " / " << endl);
177 
178  ++p0;
179  ++p1;
180  }
181  }
182 
183  STATUS("Number of differences / events: " << count[1] << " / " << count[0] << endl);
184 
185  if (buffer[0].size() != buffer[1].size()) {
186  FATAL("Different size " << buffer[0].size() << ' ' << buffer[1].size() << endl);
187  }
188 
189  if (count[1] != 0) {
190  FATAL("Number of differences " << count[1] << endl);
191  }
192 }
ROOT TTree parameter settings.
int main(int argc, char **argv)
Definition: JDiffToA.cc:117
TPaveText * p1
General purpose messaging.
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62
#define STATUS(A)
Definition: JMessage.hh:63
#define FATAL(A)
Definition: JMessage.hh:67
int debug
debug level
Definition: JSirene.cc:69
Utility class to parse command line options.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:2142
Scanning of objects from a single file according a format that follows from the extension of each fil...
Acoustic event.
Utility class to parse command line options.
Definition: JParser.hh:1698
Object reading from a list of files.
virtual bool hasNext() override
Check availability of next element.
std::ostream & operator<<(std::ostream &stream, const CLBCommonHeader &header)
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Definition: JSTDTypes.hh:14
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:448
Time-of-arrival data from acoustic piezo sensor or hydrophone.
Definition: JToA.hh:26
uint32_t SECONDS
Time of Arrival, in ns (relative to Unix epoch, 1 January 1970 00:00:00 UTC)
Definition: JToA.hh:35
uint32_t DOMID
DAQ run number.
Definition: JToA.hh:32
int64_t TOA_NS
Unique ID of the waveform that best described the signal around TOA_NS.
Definition: JToA.hh:34
uint32_t QUALITYFACTOR
The ticks (16ns) part of the DAQ frame timestamp.
Definition: JToA.hh:37
uint32_t QUALITYNORMALISATION
A measure of how good the waveform match was to the signal.
Definition: JToA.hh:38
int32_t WAVEFORMID
DOM unique identifeir.
Definition: JToA.hh:33
int32_t DETID
Definition: JToA.hh:30
int32_t RUN
detector identifier
Definition: JToA.hh:31
uint32_t TICKS
The seconds part of the DAQ frame timestamp.
Definition: JToA.hh:36
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:45