Jpp  15.0.4
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JAcousticsTriggerProcessor.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <iomanip>
3 #include <vector>
4 #include <map>
5 #include <algorithm>
6 
7 #include "TROOT.h"
8 #include "TFile.h"
9 
10 #include "JLang/JComparator.hh"
11 
12 #include "JDB/JToAshort.hh"
13 #include "JDB/JSupport.hh"
14 
15 #include "JDetector/JDetector.hh"
18 #include "JDetector/JHydrophone.hh"
19 
20 #include "JLang/JPredicate.hh"
21 #include "JLang/JComparator.hh"
22 #include "JTools/JRange.hh"
23 #include "JTools/JQuantile.hh"
24 #include "JTools/JHashMap.hh"
27 #include "JSupport/JMeta.hh"
28 
30 #include "JAcoustics/JReceiver.hh"
35 #include "JAcoustics/JEvent.hh"
37 #include "JAcoustics/JSupport.hh"
38 
39 #include "Jeep/JContainer.hh"
40 #include "Jeep/JPrint.hh"
41 #include "Jeep/JParser.hh"
42 #include "Jeep/JMessage.hh"
43 
44 #include "JTrigger/JMatch.hh"
45 #include "JTrigger/JAlgorithm.hh"
46 
47 
48 namespace JACOUSTICS {
49 
51  using JTRIGGER::JMatch;
52  using JLANG::JClonable;
53  using JTOOLS::JRange;
54 
55 
56  /**
57  * Acoustic hit.
58  */
59  struct hit_type :
60  public JPosition3D,
61  public JTransmission
62  {
63  /**
64  * Default constructor.
65  */
67  {}
68 
69 
70  /**
71  * Constructor.
72  *
73  * \param position position
74  * \param transmission transmission
75  */
76  hit_type(const JPosition3D& position, const JTransmission& transmission) :
77  JPosition3D (position),
78  JTransmission(transmission)
79  {}
80  };
81 
82 
83  /**
84  * 3D match criterion for acoustic signals.
85  */
86  class JMatch3D :
87  public JClonable< JMatch<hit_type>, JMatch3D>
88  {
89  public:
90  /**
91  * Constructor.
92  *
93  * \param V sound velocity
94  * \param Tmax_s maximal extra time [s]
95  */
96  JMatch3D(const JSoundVelocity& V, const double Tmax_s = 0.0) :
97  V(V),
98  Tmax_s(Tmax_s)
99  {}
100 
101 
102  /**
103  * Match operator.
104  *
105  * \param first hit
106  * \param second hit
107  * \return match result
108  */
109  virtual bool operator()(const hit_type& first, const hit_type& second) const override
110  {
111  using namespace JPP;
112 
113  const double t1 = V.getTime(getDistance(first.getPosition(), second.getPosition()),
114  first .getZ(),
115  second.getZ());
116 
117  const double dt = fabs(first.getToA() - second.getToA());
118 
119  return dt <= t1 + Tmax_s;
120  }
121 
122 
124  const double Tmax_s;
125  };
126 
127 
128  /**
129  * Match of two events considering overlap in time.
130  */
131  struct JEventOverlap {
132 
133  typedef JRange<double> time_range; //!< Type definition of time range.
134 
135  /**
136  * Constructor.
137  *
138  * \param Tmax_s maximal time difference between two consecutive events [s]
139  */
140  JEventOverlap(const double Tmax_s = 0.0) :
141  range(-Tmax_s, +Tmax_s)
142  {}
143 
144 
145  /**
146  * Match criterion.
147  *
148  * \param first first event
149  * \param second second event
150  * \return true if two events overlap in time; else false
151  */
152  bool operator()(const JEvent& first, const JEvent& second) const
153  {
154  using namespace JPP;
155 
156  if (first .empty()) return false;
157  if (second.empty()) return false;
158 
159  time_range r1(make_array(first .begin(), first .end(), &JTransmission::getToA));
160  time_range r2(make_array(second.begin(), second.end(), &JTransmission::getToA));
161 
162  r1.add(range);
163  r2.add(range);
164 
165  return overlap(r1, r2);
166  }
167 
169  };
170 }
171 
172 
173 /**
174  * \file
175  *
176  * Main program to trigger acoustic data.
177  *
178  * An acoustic event is based on coincidences between times of arrival.\n
179  * If the number of coincident times of arrival exceeds a preset minimum,
180  * the event is triggered and subsequently stored in the output file.
181  * \author mdejong
182  */
183 int main(int argc, char **argv)
184 {
185  using namespace std;
186  using namespace JPP;
187 
189 
190  typedef JContainer< vector<JHydrophone> > hydrophones_container;
191 
193  JLimit_t& numberOfEvents = inputFile.getLimit();
196  JSoundVelocity V = getSoundVelocity; // default sound velocity
197  string detectorFile;
198  hydrophones_container hydrophones; // hydrophones
199  double precision;
200  int ID;
201  set<int> waveforms;
202  int debug;
203 
204  try {
205 
206  JParser<> zap("Main program to trigger acoustic data.");
207 
208  zap['f'] = make_field(inputFile, "output of JConvertDB -q toashort");
209  zap['n'] = make_field(numberOfEvents) = JLimit::max();
210  zap['@'] = make_field(parameters, "trigger parameters");
211  zap['o'] = make_field(outputFile, "output file") = "event.root";
212  zap['V'] = make_field(V, "sound velocity") = JPARSER::initialised();
213  zap['a'] = make_field(detectorFile, "detector file");
214  zap['H'] = make_field(hydrophones, "hydrophone data") = JPARSER::initialised();
215  zap['p'] = make_field(precision, "precision time-of-arrival") = 1.0e-6;
216  zap['D'] = make_field(ID) = 10000;
217  zap['w'] = make_field(waveforms) = JPARSER::initialised();
218  zap['d'] = make_field(debug) = 1;
219 
220  zap(argc, argv);
221  }
222  catch(const exception &error) {
223  FATAL(error.what() << endl);
224  }
225 
226 
228 
229  try {
230  load(detectorFile, detector);
231  }
232  catch(const JException& error) {
233  FATAL(error);
234  }
235 
236  V.set(detector.getUTMZ());
237 
238 
239  const int FACTORY_LIMIT = 10000;
240  const double TMaxExtra_s = 1.0e-3;
241 
242  const JMatch3D match(V, TMaxExtra_s);
243  const JEventOverlap overlap; // (parameters.TMax_s);
244 
245 
246  JHashMap<int, JReceiver> receivers;
247 
248  for (JDetector::const_iterator i = detector.begin(); i != detector.end(); ++i) {
249 
250  JPosition3D pos(0.0, 0.0, 0.0);
251 
252  if (i->getFloor() == 0) { // get relative position of hydrophone
253 
254  try {
255  pos = getPosition(hydrophones.begin(),
256  hydrophones.end(),
257  make_predicate(&JHydrophone::getString, i->getString()));
258  }
259  catch(const exception&) {
260  continue; // if no position avialable, discard hydrophone
261  }
262  }
263 
264  receivers[i->getID()] = JReceiver(i->getID(),
265  i->getPosition() + pos,
266  i->getT0() * 1.0e-9);
267  }
268 
269 
270  outputFile.open();
271 
272  if (!outputFile.is_open()) {
273  FATAL("Error opening file " << outputFile << endl);
274  }
275 
276  outputFile.put(JMeta(argc, argv));
277  outputFile.put(parameters);
278 
279 
280  // input data
281 
282  string oid; // detector identifier
283 
284  typedef vector<hit_type> buffer_type; // data type
285 
286  map<int, buffer_type> f1; // receiver -> data
287 
288  while (inputFile.hasNext()) {
289 
290  if (inputFile.getCounter()%1000 == 0) {
291  STATUS("counter: " << setw(8) << inputFile.getCounter() << '\r' << flush); DEBUG(endl);
292  }
293 
294  JToAshort* parameters = inputFile.next();
295 
296  if (oid == "") {
297  oid = parameters->DETID;
298  }
299 
300  if (oid != parameters->DETID) { // consistency check
301  FATAL("Invalid detector identifier " << parameters->DETID << " != " << oid << endl);
302  }
303 
304  if (waveforms.empty() || waveforms.count(parameters->EMITTERID)) {
305 
306  if (receivers.has(parameters->DOMID)) {
307 
308  const JReceiver& receiver = receivers[parameters->DOMID];
309 
310  const JTransmission transmission(parameters->RUNNUMBER,
311  parameters->DOMID,
312  parameters->QUALITYFACTOR,
313  parameters->UNIXTIMEBASE + receiver.getT(parameters->TOA_S),
314  parameters->UNIXTIMEBASE + receiver.getT(parameters->TOA_S));
315 
316  f1[receiver.getID()].push_back(hit_type(receiver.getPosition(), transmission));
317  }
318  }
319  }
320  STATUS(endl);
321 
322 
323  // filter similar hits
324 
325  for (map< int, buffer_type>::iterator receiver = f1.begin(); receiver != f1.end(); ++receiver) {
326 
327  buffer_type& buffer = receiver->second;
328 
329  sort(buffer.begin(), buffer.end(), JTransmission::compare(precision));
330 
331  buffer.erase(unique(buffer.begin(), buffer.end(), JTransmission::equals(precision)), buffer.end());
332  }
333 
334 
335  buffer_type data;
336 
337  for (map<int, buffer_type>::iterator receiver = f1.begin(); receiver != f1.end(); ++receiver) {
338 
339  // selection based on quality
340 
341  buffer_type& f2 = receiver->second;
342  double Qmin = parameters.Q;
343 
344  if (parameters.Q > 0.0 && parameters.Q < 1.0) {
345 
346  const JQuantile Q1("quality", make_array(f2.begin(), f2.end(), &JTransmission::getQ), true);
347 
348  Qmin = Q1.getQuantile(parameters.Q);
349  }
350 
351  buffer_type::iterator __end = partition(f2.begin(), f2.end(), make_predicate(&JTransmission::getQ, Qmin, JComparison::gt()));
352 
353  copy(f2.begin(), __end, back_inserter(data));
354  }
355 
356 
357  sort(data.begin(), data.end(), make_comparator(&hit_type::getToA, JComparison::lt())); // sort according time-of-arrival
358 
359 
360  buffer_type buffer(FACTORY_LIMIT); // local buffer of hits
361  JEvent out[2]; // FIFO of events
362  int number_of_events = 0;
363 
364  for (buffer_type::const_iterator p = data.begin(); p != data.end(); ++p) {
365 
366  if (distance(data.cbegin(),p)%1000 == 0) {
367  STATUS("processed: " << FIXED(5,1) << (double) (100 * distance(data.cbegin(),p)) / (double) data.size() << "%" << '\r' << flush); DEBUG(endl);
368  }
369 
370  buffer_type::const_iterator q = p;
371 
372  while (++q != data.end() && q->getToA() - p->getToA() <= parameters.TMax_s) {}
373 
374  if (distance(p,q) >= parameters.numberOfHits) {
375 
376  if (distance(p,q) < FACTORY_LIMIT) {
377 
378  buffer_type::iterator root = buffer.begin();
379  buffer_type::iterator __p = buffer.begin();
380  buffer_type::iterator __q = buffer.begin();
381 
382  *root = *p;
383 
384  ++__p;
385  ++__q;
386 
387  for (buffer_type::const_iterator i = p; ++i != q; ) {
388  if (match(*p,*i)) {
389  *__q = *i;
390  ++__q;
391  }
392  }
393 
394  if (distance(root,__q) >= parameters.numberOfHits) {
395 
396  __q = clusterize(__p, __q, match, parameters.numberOfHits - 1);
397 
398  if (distance(root,__q) >= parameters.numberOfHits) {
399  out[1] = JEvent(oid, ++number_of_events, ID, root, __q);
400  }
401  }
402 
403  } else {
404 
405  out[1] = JEvent(oid, ++number_of_events, ID, p, q);
406  }
407 
408  if (!overlap(out[0],out[1])) {
409 
410  if (!out[0].empty()) {
411  outputFile.put(out[0]); // write
412  }
413 
414  out[0] = out[1]; // shift
415 
416  } else {
417 
418  out[0].merge(out[1]); // merge
419  }
420 
421  out[1].clear();
422  }
423  }
424 
425  if (!out[0].empty()) {
426  outputFile.put(out[0]); // write
427  }
428  STATUS(endl);
429 
430  STATUS("triggers: " << setw(7) << number_of_events << endl);
431 
432  JMultipleFileScanner<JMeta> io(inputFile);
433 
434  io >> outputFile;
435 
436  outputFile.close();
437 }
static struct JTRIGGER::@74 clusterize
Anonymous structure for clustering of hits.
Auxiliary class for ROOT I/O of application specific meta data.
Definition: JMeta.hh:70
hit_type(const JPosition3D &position, const JTransmission &transmission)
Constructor.
Object writing to file.
double UNIXTIMEBASE
[s]
Definition: JToAshort.hh:25
Utility class to parse command line options.
Definition: JParser.hh:1500
General exception.
Definition: JException.hh:23
Acoustic receiver.
Definition: JReceiver.hh:27
JPredicate< JResult_t T::*, JComparison::eq > make_predicate(JResult_t T::*member, const JResult_t value)
Helper method to create predicate for data member.
Definition: JPredicate.hh:128
int main(int argc, char *argv[])
Definition: Main.cc:15
double getQ() const
Get quality.
Greater than.
Definition: JComparison.hh:73
JComparator< JResult_t T::*, JComparison::lt > make_comparator(JResult_t T::*member)
Helper method to create comparator between values of data member.
Definition: JComparator.hh:185
Sound velocity.
Algorithms for hit clustering and sorting.
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
JVector3D getPosition(T __begin, T __end, const JPredicate< JTypename_t, JComparator_t > &predicate)
Get position from element in data which corresponds to given predicate.
General purpose class for hash map of unique keys.
Definition: JHashMap.hh:72
Auxiliary data structure for running average, standard deviation and quantiles.
Definition: JQuantile.hh:43
General purpose class for hash map of unique elements.
#define STATUS(A)
Definition: JMessage.hh:63
ROOT TTree parameter settings.
Detector data structure.
Definition: JDetector.hh:89
then JPlot1D f $WORKDIR postfit[prefit\] root
Recording of objects on file according a format that follows from the file name extension.
virtual bool operator()(const hit_type &first, const hit_type &second) const override
Match operator.
std::string DETID
constraint
Definition: JToAshort.hh:22
*fatal Wrong number of arguments esac JCookie sh typeset Z DETECTOR typeset Z SOURCE_RUN typeset Z TARGET_RUN set_variable PARAMETERS_FILE $WORKDIR parameters
Definition: diff-Tuna.sh:38
Acoustic event.
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition: JParser.hh:66
then let ID
Definition: JAcoustics.sh:31
Function object interface for hit matching.
Definition: JMatch.hh:25
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:446
3D match criterion for acoustic signals.
static const int FACTORY_LIMIT
Bit indicating max nhits reached in trigger.
Definition: trigger.hh:16
JMatch3D(const JSoundVelocity &V, const double Tmax_s=0.0)
Constructor.
string outputFile
Data structure for detector geometry and calibration.
then echo The file $DIR KM3NeT_00000001_00000000 root already please rename or remove it first
double getQuantile(const double Q, const bool reverse=false) const
Get quantile.
Definition: JQuantile.hh:309
double getDistance(const JFirst_t &first, const JSecond_t &second)
Get distance between objects.
Data structure for hydrophone.
bool operator()(const JEvent &first, const JEvent &second) const
Match criterion.
Type list.
Definition: JTypeList.hh:22
JFIT::JEvent JEvent
Definition: JHistory.hh:282
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
static const JSoundVelocity getSoundVelocity(1541.0,-17.0e-3,-2000.0)
Function object for velocity of sound.
Base class for match operations for cluster and hit-preprocessing methods.
Acoustic transmission.
I/O formatting auxiliaries.
Detector file.
Definition: JHead.hh:196
Acoustics support kit.
Detector support kit.
Acoustics toolkit.
Auxiliary wrapper for I/O of container with optional comment (see JComment).
Definition: JContainer.hh:39
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1961
range_type & add(argument_type x)
Add offset.
Definition: JRange.hh:447
const array_type< JValue_t > & make_array(const JValue_t(&array)[N])
Method to create array of values.
Definition: JVectorize.hh:54
int getID() const
Get identifier.
Definition: JObjectID.hh:50
ROOT I/O of application specific meta data.
Template class for object cloning.
Definition: JClonable.hh:20
ROOT TTree parameter settings.
int debug
debug level
Definition: JSirene.cc:63
const JPosition3D & getPosition() const
Get position.
Definition: JPosition3D.hh:130
Range of values.
Definition: JRange.hh:38
General purpose messaging.
Implementation for velocity of sound.
#define FATAL(A)
Definition: JMessage.hh:67
Scanning of objects from multiple files according a format that follows from the extension of each fi...
double getT(const double t_s) const
Get corrected time.
Definition: JReceiver.hh:72
Acoustic receiver.
z range($ZMAX-$ZMIN)< $MINIMAL_DZ." fi fi typeset -Z 4 STRING typeset -Z 2 FLOOR JPlot1D -f $
JEventOverlap(const double Tmax_s=0.0)
Constructor.
virtual double getTime(const double D_m, const double z1, const double z2) const override
Get propagation time of sound.
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
Acoustic trigger parameters.
Auxiliary class to define a range between two values.
General purpose class for object reading from a list of file names.
Utility class to parse command line options.
hit_type()
Default constructor.
double getToA() const
Get calibrated time of arrival.
then usage $script< input_file >< detector_file > fi set_variable OUTPUT_DIR set_variable SELECTOR JDAQTimesliceL1 set_variable DEBUG case set_variable DEBUG
Acoustic transmission.
Auxiliary class to compare transmissions.
void copy(const Head &from, JHead &to)
Copy header from from to to.
Definition: JHead.cc:139
Acoustic event.
double TOA_S
[s]
Definition: JToAshort.hh:28
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:73
Data structure for position in three dimensions.
Definition: JPosition3D.hh:36
bool equals(const JFirst_t &first, const JSecond_t &second, const double precision=std::numeric_limits< double >::min())
Check equality.
Definition: JMathToolkit.hh:86
bool overlap(const JRange< T, JComparator_t > &first, const JRange< T, JComparator_t > &second)
Test overlap between ranges.
Definition: JRange.hh:663
JRange< double > time_range
Type definition of time range.
do set_variable DETECTOR_TXT $WORKDIR detector
JSoundVelocity & set(const double z0)
Set depth.
Match of two events considering overlap in time.
double getZ() const
Get z position.
Definition: JVector3D.hh:115
int EMITTERID
waveform identifier
Definition: JToAshort.hh:27
3D match criterion.
Definition: JMatch3D.hh:29
Container I/O.