Jpp  master_rocky-43-ge265d140c
the software that should make you happy
JAcousticsEventMonitor.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <iomanip>
3 
4 #include "JAcoustics/JEvent.hh"
5 #include "JAcoustics/JSupport.hh"
6 
7 #include "JTools/JQuantile.hh"
8 
10 #include "JSupport/JTreeScanner.hh"
11 
12 #include "Jeep/JPrint.hh"
13 #include "Jeep/JParser.hh"
14 #include "Jeep/JMessage.hh"
15 
16 namespace {
17 
19 
20  /**
21  * Auxiliary data structure to compare transmissions.
22  */
23  struct JCompare {
24  /**
25  * Compare transmissions.
26  *
27  * \param first first transmission
28  * \param second second transmission
29  * \return true if first transmission before second; else false
30  */
31  inline bool operator()(const JTransmission& first, const JTransmission& second) const
32  {
33  if (first.getID() == second.getID())
34  if (fabs(first.getToA() - second.getToA()) <= waveform_s)
35  return false;
36  else
37  return first.getToA() < second.getToA();
38  else
39  return first.getID() < second.getID();
40  }
41 
42  double waveform_s;
43 
44  } compare;
45 }
46 
47 
48 /**
49  * \file
50  *
51  * Example program to monitor acoustic events.
52  * \author mdejong
53  */
54 int main(int argc, char **argv)
55 {
56  using namespace std;
57  using namespace JPP;
58 
59  JMultipleFileScanner<> inputFile;
60  JLimit_t& numberOfEvents = inputFile.getLimit();
61  double Tmax_s;
62  int debug;
63 
64  try {
65 
66  JParser<> zap("Example program to monitor acoustic events.");
67 
68  zap['f'] = make_field(inputFile, "output of JAcousticEventBuilder[.sh]");
69  zap['n'] = make_field(numberOfEvents) = JLimit::max();
70  zap['p'] = make_field(compare.waveform_s, "wavefrom duration [s]") = 5.0e-3;
71  zap['D'] = make_field(Tmax_s, "time difference between events [s]") = 20.0e-3;
72  zap['d'] = make_field(debug) = 2;
73 
74  zap(argc, argv);
75  }
76  catch(const exception &error) {
77  FATAL(error.what() << endl);
78  }
79 
80 
82 
83  JEvent* p0 = new JEvent();
84 
85  counter_type events = 0;
86 
87  JQuantile Q;
88 
89  for (JMultipleFileScanner<>::const_iterator i = inputFile.begin(); i != inputFile.end(); ++i) {
90 
91  JTreeScanner_t in(*i);
92 
93  for (JTreeScanner_t::iterator p1 = in.begin(); p1 != in.end() && events != inputFile.getLimit(); ++p1, ++events) {
94 
95  STATUS("event " << setw(8) << events << '\r'); DEBUG(endl);
96 
97  if (!p1->empty() && !p0->empty()) {
98 
99  if (p1->begin()->getToE() < p0->rbegin()->getToE() + Tmax_s) {
100 
101  if (debug >= debug_t) {
102 
103  cout << endl;
104 
105  cout << p0->getDetectorID() << endl
106  << FIXED(20,6) << p0-> begin()->getToE() << ' '
107  << FIXED(20,6) << p0->rbegin()->getToE() << endl
108  << FIXED(20,6) << p1-> begin()->getToE() << ' '
109  << FIXED(20,6) << p1->rbegin()->getToE() << endl;
110  }
111 
112  JEvent evt[] = { *p0, *p1 };
113 
114  for (int i = 0; i != 2; ++i) {
115  sort(evt[i].begin(), evt[i].end(), compare);
116  }
117 
118  size_t n = 0;
119 
120  for (JEvent::const_iterator
121  i0 = evt[0].begin(),
122  i1 = evt[1].begin(); i0 != evt[0].end() && i1 != evt[1].end(); ) {
123 
124  if (compare(*i0, *i1))
125  ++i0;
126  else if (compare(*i1, *i0))
127  ++i1;
128 
129  else {
130 
131  n += 1;
132 
133  if (debug >= debug_t) {
134 
135  cout << setw(2) << p0->getID() << ' '
136  << setw(10) << i0->getID() << ' '
137  << FIXED(20,6) << i0->getToA() << ' '
138  << FIXED( 8,0) << i0->getQ() << ' ';
139 
140  cout << setw(2) << p1->getID() << ' '
141  << setw(10) << i1->getID() << ' '
142  << FIXED(20,6) << i1->getToA() << ' '
143  << FIXED( 8,0) << i1->getQ() << ' ';
144 
145  cout << FIXED( 9,6) << (i1->getToA() - i0->getToA()) << endl;
146  }
147 
148  if (i0->getToA() < i1->getToA())
149  ++i0;
150  else if (i1->getToA() < i0->getToA())
151  ++i1;
152  else {
153  ++i0;
154  ++i1;
155  }
156  }
157  }
158 
159  if (n != 0) {
160  Q.put(n);
161  }
162  }
163  }
164 
165  new (p0) JEvent(*p1);
166  }
167  }
168  STATUS(endl);
169 
170  delete p0;
171 
172  cout << "Number of errors / events: " << Q.getCount() << " / " << events << endl;
173 
174  Q.print(cout);
175 }
int main(int argc, char **argv)
Acoustic event.
ROOT TTree parameter settings.
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
Scanning of objects from multiple files according a format that follows from the extension of each fi...
Utility class to parse command line options.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:2142
I/O formatting auxiliaries.
Utility class to parse command line options.
Definition: JParser.hh:1698
Base class for JTreeScanner.
Definition: JTreeScanner.hh:56
Template definition for direct access of elements in ROOT TChain.
@ debug_t
debug
Definition: JMessage.hh:29
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
JFIT::JEvent JEvent
Definition: JHistory.hh:353
Long64_t counter_type
Type definition for counter.
const int n
Definition: JPolint.hh:786
Definition: JSTDTypes.hh:14
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:448
int getID() const
Get emitter identifier.
const int getDetectorID() const
Get detector identifier.
Acoustic transmission.
double getToA() const
Get calibrated time of arrival.
int getID() const
Get identifier.
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:45
Auxiliary data structure for running average, standard deviation and quantiles.
Definition: JQuantile.hh:46
std::ostream & print(std::ostream &out, bool lpr=true) const
Print quantile.
Definition: JQuantile.hh:382
void put(const double x, const double w=1.0)
Put value.
Definition: JQuantile.hh:133
long long int getCount() const
Get total count.
Definition: JQuantile.hh:186