Jpp 20.0.0-rc.9-29-gccc23c492-D
the software that should make you happy
Loading...
Searching...
No Matches
JAcousticsTriggerProcessor.cc
Go to the documentation of this file.
1#include <string>
2#include <iostream>
3#include <iomanip>
4#include <vector>
5#include <set>
6#include <map>
7#include <algorithm>
8
9#include "TROOT.h"
10#include "TFile.h"
11#include "TH1D.h"
12#include "TGraph.h"
13
14#include "JDetector/JDetector.hh"
15#include "JDetector/JDetectorToolkit.hh"
16#include "JDetector/JHydrophone.hh"
17
18#include "JLang/JClonable.hh"
19#include "JLang/JPredicate.hh"
20#include "JLang/JComparator.hh"
21
22#include "JTools/JHashMap.hh"
23
24#include "JROOT/JManager.hh"
25#include "JROOT/JGraph.hh"
26
27#include "JSupport/JMultipleFileScanner.hh"
28#include "JSupport/JFileRecorder.hh"
29#include "JSupport/JMeta.hh"
30
31#include "JAcoustics/JToA.hh"
39#include "JAcoustics/JEvent.hh"
41
42#include "JGeometry3D/JVector3D.hh"
43#include "JGeometry3D/JVertex3D.hh"
44#include "JGeometry3D/JPosition3D.hh"
45#include "JGeometry3D/JCylinder3D.hh"
46
47#include "Jeep/JContainer.hh"
48#include "Jeep/JProperties.hh"
49#include "Jeep/JPrint.hh"
50#include "Jeep/JParser.hh"
51#include "Jeep/JMessage.hh"
52
53#include "JTrigger/JMatch.hh"
54#include "JTrigger/JAlgorithm.hh"
55
56
57namespace JACOUSTICS {
58
59 using JGEOMETRY3D::JVector3D;
60 using JGEOMETRY3D::JVertex3D;
61 using JGEOMETRY3D::JPosition3D;
62 using JTRIGGER::JMatch;
63 using JLANG::JClonable;
64
65
66 /**
67 * Transmission with position.
68 */
69 struct hit_type :
70 public JPosition3D,
71 public JTransmission
72 {
73 /**
74 * Default constructor.
75 */
77 {}
78
79
80 /**
81 * Constructor.
82 *
83 * \param position position
84 * \param transmission transmission
85 */
86 hit_type(const JPosition3D& position, const JTransmission& transmission) :
87 JPosition3D (position),
88 JTransmission(transmission)
89 {}
90 };
91
92
93 /**
94 * 3D match criterion for acoustic signals.
95 */
96 class JMatch3D :
97 public JClonable< JMatch<hit_type>, JMatch3D>
98 {
99 public:
100 /**
101 * Constructor.
102 *
103 * \param V sound velocity
104 * \param Tmax_s maximal extra time [s]
105 */
106 JMatch3D(const JAbstractSoundVelocity& V, const double Tmax_s = 0.0) :
107 V(V),
109 {}
110
111
112 /**
113 * Match operator.
114 *
115 * \param first hit
116 * \param second hit
117 * \return match result
118 */
119 virtual bool operator()(const hit_type& first, const hit_type& second) const override
120 {
121 using namespace JPP;
122
123 const double t1 = V.getTime(getDistance(first.getPosition(), second.getPosition()),
124 first .getZ(),
125 second.getZ());
126 const double dt = fabs(first.getToA() - second.getToA());
127
128 return dt <= t1 + Tmax_s;
129 }
130
131
133 const double Tmax_s;
134 };
135
136
137 /**
138 * Event with vertex.
139 */
140 struct event_type :
141 public JVertex3D,
142 public JEvent
143 {
144 /**
145 * Default constructor.
146 */
148 JVertex3D(),
149 JEvent()
150 {}
151
152
153 /**
154 * Copy constructor.
155 *
156 * \param event event
157 */
158 event_type(const JEvent& event) :
159 JVertex3D(),
160 JEvent(event)
161 {}
162
163
164 /**
165 * Constructor.
166 *
167 * \param vertex vertex
168 * \param event event
169 */
170 event_type(const JVertex3D& vertex,
171 const JEvent& event) :
172 JVertex3D(vertex),
173 JEvent(event)
174 {}
175
176
177 /**
178 * Get average quality.
179 *
180 * \return quality
181 */
182 inline double getQ() const
183 {
184 double Q = 0.0;
185
186 if (!this->empty()) {
187
188 for (const_iterator i = this->begin(); i != this->end(); ++i) {
189 Q += i->getQ();
190 }
191
192 Q /= this->size();
193 }
194
195 return Q;
196 }
197
198
199 /**
200 * Write event to output stream.
201 *
202 * \param out output stream
203 * \return object event
204 * \return output stream
205 */
206 friend inline std::ostream& operator<<(std::ostream& out, const event_type& object)
207 {
208 using namespace std;
209
210 return out << setw(3) << object.getID() << ' '
211 << FIXED(9,3) << object.getT() - TOA_s << " [s]" << ' '
212 << setw(4) << object.size() << ' '
213 << FIXED(9,2) << object.getQ() << ' '
214 << object.getPosition();
215 }
216
217 static double TOA_s; //!< start time of data
218 };
219
220 /**
221 * Set start time of data.
222 */
223 double event_type::TOA_s = std::numeric_limits<double>::max();
224
225
226 /**
227 * Match of two events considering overlap in time and position.
228 */
230 /**
231 * Constructor.
232 *
233 * \param Tmax_s maximal time difference between two events [s]
234 * \param Dmax_m maximal distance between two events [m]
235 */
236 JEventOverlap(const double Tmax_s, const double Dmax_m = std::numeric_limits<double>::max()) :
237 Tmax_s(Tmax_s),
239 {}
240
241
242 /**
243 * Match criterion.
244 *
245 * \param first first event
246 * \param second second event
247 * \return true if two events overlap; else false
248 */
249 bool operator()(const event_type& first, const event_type& second) const
250 {
251 using namespace JPP;
252
253 return fabs(first.getT() - second.getT()) <= Tmax_s && getDistance(first.getPosition(), second.getPosition()) <= Dmax_m;
254 }
255
256 const double Tmax_s;
257 const double Dmax_m;
258 };
259
260
261 /**
262 * Vertex.
263 */
264 struct vertex_type :
265 public JVertex3D
266 {
267 /**
268 * Default constructor.
269 */
271 JVertex3D(),
272 N(0),
273 Q(0.0)
274 {}
275
276
277 /**
278 * Constructor.
279 *
280 * \param p0 position
281 * \param t0 time
282 * \param n0 number of hits
283 * \param q0 total quality
284 */
285 vertex_type(const JVector3D& p0, const double t0, const int n0, const double q0) :
286 JVertex3D(p0, t0),
287 N(n0),
288 Q(q0)
289 {}
290
291
292 int N; //!< number of hits
293 double Q; //!< total quality
294 };
295
296
297 /**
298 * Vertex locator.
299 */
300 struct JVelo :
301 public std::vector<JPosition3D>
302 {
303 /**
304 * Constructor.
305 *
306 * \param V sound velocity
307 * \param RMax_m radius for generation of vertices [m]
308 * \param Xv_m step size for generation of vertices [m]
309 * \param factor multiplication factor for corresponding time window
310 */
312 const double RMax_m,
313 const double Xv_m,
314 const double factor = 1.0) :
315 V(V),
316 Tmax_s(factor * Xv_m / V())
317 {
318 this->push_back(JVector3D(0.0,0.0,0.0));
319
320 if (RMax_m > 0.0 && Xv_m > 0.0) {
321 for (double x = 0.5*Xv_m; x <= RMax_m; x += Xv_m) {
322 for (double y = 0.5*Xv_m; y <= RMax_m; y += Xv_m) {
323 if (x*x + y*y <= RMax_m*RMax_m) {
324 this->push_back(JVector3D(-x, +y, 0.0));
325 this->push_back(JVector3D(+x, +y, 0.0));
326 this->push_back(JVector3D(-x, -y, 0.0));
327 this->push_back(JVector3D(+x, -y, 0.0));
328 }
329 }
330 }
331 }
332 }
333
334
335 /**
336 * Check vertex.
337 *
338 * \param vx vertex
339 * \param __begin begin of data
340 * \param __end end of data
341 * \return number of hits
342 */
343 template<class T>
344 inline int operator()(const JVertex3D& vx, T __begin, T __end) const
345 {
346 const JPosition3D& p0 = vx.getPosition();
347 const double t0 = vx.getT();
348
349 int n0 = 0;
350
351 for (T p = __begin; p != __end; ++p) {
352
353 const double t1 = p->getToA() - V.getTime(p0.getDistance(p->getPosition()), p0.getZ(), p->getZ());
354
355 if (fabs(t1 - t0) <= Tmax_s) {
356 n0 += 1;
357 }
358 }
359
360 return n0;
361 }
362
363
364 /**
365 * Locate vertex around given position.
366 *
367 * \param position position
368 * \param root root hit
369 * \param __begin begin of data
370 * \param __end end of data
371 * \param numberOfHits number of hits
372 * \return vertex
373 */
374 template<class T>
375 inline vertex_type operator()(const JPosition3D& position, const hit_type& root, T __begin, T __end, const int numberOfHits = 0) const
376 {
377 vertex_type vertex;
378
379 for (const_iterator i = this->cbegin(); i != this->cend(); ++i) {
380
381 const JPosition3D p0 = position + *i;
382 const double d0 = p0.getDistance(root.getPosition());
383 const double t0 = root.getToA() - V.getTime(d0, p0.getZ(), root.getZ());
384
385 int n0 = 1;
386 double q0 = root.getQ();
387
388 for (T p = __begin; p != __end && n0 + distance(p, __end) >= numberOfHits && n0 + distance(p, __end) >= vertex.N; ++p) {
389
390 const double d1 = p0.getDistance(p->getPosition());
391 const double t1 = p->getToA() - V.getTime(d1, p0.getZ(), p->getZ());
392
393 if (fabs(t1 - t0) <= Tmax_s) {
394 n0 += 1;
395 q0 += p->getQ();
396 }
397 }
398
399 if (q0 > vertex.Q) {
400 vertex = vertex_type(p0, t0, n0, q0);
401 }
402 }
403
404 return vertex;
405 }
406
408 const double Tmax_s;
409 };
410}
411
412
413/**
414 * \file
415 *
416 * Main program to trigger acoustic data.
417 *
418 * An acoustic event is based on coincidences between times of arrival.\n
419 * If the number of coincident times of arrival exceeds a preset minimum,
420 * the event is triggered and subsequently stored in the output file.
421 * \author mdejong
422 */
423int main(int argc, char **argv)
424{
425 using namespace std;
426 using namespace JPP;
427
428 typedef JTYPELIST<JEvent, JTriggerParameters, JMeta, TH1D, TGraph>::typelist typelist;
429
430 typedef JContainer< vector<JHydrophone> > hydrophones_container;
431
432 JMultipleFileScanner<JToA> inputFile;
433 JLimit_t& numberOfEvents = inputFile.getLimit();
434 JTriggerParameters parameters;
435 int factoryLimit = 10000;
436 double TMaxExtra_s = 500.0e-6;
437 double RMax_m = 0.0;
438 double DMax_m = 0.0;
439 double Xv_m = 0.0;
440 double factor = 2.0;
441 double Z_m = 0.0; // z-position of emitters
442 double fraction = 0.75; // fraction of working modules
443 double grid = 10.0;
444 JFileRecorder<typelist> outputFile;
445 JSoundVelocity V = getSoundVelocity; // default sound velocity
446 string detectorFile;
447 hydrophones_container hydrophones; // hydrophones
448 double precision;
449 set<int> ID;
450 int debug;
451
452 try {
453
454 JProperties properties;
455
456 properties.insert(gmake_property(parameters.Q));
457 properties.insert(gmake_property(parameters.TMax_s));
458 properties.insert(gmake_property(parameters.numberOfHits));
459 properties.insert(gmake_property(fraction));
460 properties.insert(gmake_property(factoryLimit));
461 properties.insert(gmake_property(TMaxExtra_s));
462 properties.insert(gmake_property(RMax_m));
463 properties.insert(gmake_property(DMax_m));
464 properties.insert(gmake_property(Xv_m));
465 properties.insert(gmake_property(factor));
466 properties.insert(gmake_property(Z_m));
467 properties.insert(gmake_property(grid));
468
469 JParser<> zap("Main program to trigger acoustic data.");
470
471 zap['f'] = make_field(inputFile, "output of JToA");
472 zap['n'] = make_field(numberOfEvents) = JLimit::max();
473 zap['@'] = make_field(properties, "trigger parameters");
474 zap['o'] = make_field(outputFile, "output file") = "event.root";
475 zap['V'] = make_field(V, "sound velocity") = JPARSER::initialised();
476 zap['a'] = make_field(detectorFile, "detector file");
477 zap['H'] = make_field(hydrophones, "hydrophone data") = JPARSER::initialised();
478 zap['p'] = make_field(precision, "precision time-of-arrival") = 1.0e-6;
479 zap['W'] = make_field(ID, "waveform identifiers") = JPARSER::initialised();
480 zap['d'] = make_field(debug) = 1;
481
482 zap(argc, argv);
483 }
484 catch(const exception &error) {
485 FATAL(error.what() << endl);
486 }
487
488
489 JDetector detector;
490
491 try {
492 load(detectorFile, detector);
493 }
494 catch(const JException& error) {
495 FATAL(error);
496 }
497
498 V.set(detector.getUTMZ());
499
500 const JCylinder3D cylinder(detector.begin(), detector.end());
501
502 const JVector3D center(cylinder.getX(), cylinder.getY(), Z_m);
503
504 //const JConstantSoundVelocity V0 = V(cylinder.getZmax());
505 const JAbstractSoundVelocity& V0 = V;
506
507 const JVelo Velo(V0, RMax_m, Xv_m*grid, factor);
508 const JVelo velo(V0, RMax_m/grid, Xv_m, factor);
509
510 const JMatch3D match(V0, TMaxExtra_s);
511 const JEventOverlap overlap2D(parameters.TMax_s, DMax_m);
512 const JEventOverlap overlap1D(parameters.TMax_s);
513
514 JHashMap<int, JReceiver> receivers;
515
516 for (JDetector::const_iterator i = detector.begin(); i != detector.end(); ++i) {
517
518 JPosition3D pos(0.0, 0.0, 0.0);
519
520 if (i->getFloor() == 0) { // get relative position of hydrophone
521
522 try {
523 pos = getPosition(hydrophones.begin(),
524 hydrophones.end(),
525 make_predicate(&JHydrophone::getString, i->getString()));
526 }
527 catch(const exception&) {
528 continue; // if no position available, discard hydrophone
529 }
530 }
531
532 receivers[i->getID()] = JReceiver(i->getID(),
533 i->getPosition() + pos,
534 i->getT0() * 1.0e-9);
535 }
536
537 // monitoring
538
539 JManager<int, TH1D> M1(new TH1D("M1[%]", NULL, 2001, -0.5, 2000.5));
540 JManager<int, TH1D> M2(new TH1D("M2[%]", NULL, 2001, -0.5, 2000.5));
541 JManager<int, TH1D> M3(new TH1D("M3[%]", NULL, 2001, -0.5, 2000.5));
542 JManager<int, TH1D> M4(new TH1D("M4[%]", NULL, 2001, -0.5, 2000.5));
543 JManager<int, TH1D> M5(new TH1D("M5[%]", NULL, 2001, -0.5, 2000.5));
544
545 map<int, JGraph_t> G2;
546
547 outputFile.open();
548
549 if (!outputFile.is_open()) {
550 FATAL("Error opening file " << outputFile << endl);
551 }
552
553 outputFile.put(JMeta(argc, argv));
554 outputFile.put(parameters);
555
556
557 // input data
558
559 typedef vector<hit_type> buffer_type; // data type
560
561 map<int, map<int, buffer_type> > f1; // emitter -> receiver -> data
562
563 set<int> ids; // working receivers
564
565 STATUS("reading data... " << flush); DEBUG(endl);
566
567 while (inputFile.hasNext()) {
568
569 JToA* parameters = inputFile.next();
570
571 if (detector.getID() != parameters->DETID) { // consistency check
572 FATAL("Invalid detector identifier " << parameters->DETID << " != " << detector.getID() << endl);
573 }
574
575 if (receivers.has(parameters->DOMID)) {
576
577 ids.insert(parameters->DOMID);
578
579 const JReceiver& receiver = receivers[parameters->DOMID];
580 const double toa_s = parameters->TOA_S();
581
582 if (toa_s < event_type::TOA_s) {
583 event_type::TOA_s = toa_s;
584 }
585
586 const JTransmission transmission(parameters->RUN,
587 parameters->DOMID,
588 parameters->QUALITYFACTOR,
589 parameters->QUALITYNORMALISATION,
590 receiver.getT(toa_s),
591 receiver.getT(toa_s));
592
593 if (ID.empty() || ID.count(getWaveformID(parameters->WAVEFORMID)) != 0) {
594 f1[getWaveformID(parameters->WAVEFORMID)][receiver.getID()].push_back(hit_type(receiver.getPosition(), transmission));
595 }
596 }
597 }
598 STATUS("OK" << endl);
599
600 if (parameters.numberOfHits == 0) {
601
602 parameters.numberOfHits = (int) (ids.size() * fraction);
603
604 STATUS("Number of hits " << parameters.numberOfHits << endl);
605 }
606
607 vector<event_type> queue;
608
609 for (map<int, map<int, buffer_type> >::iterator i = f1.begin(); i != f1.end(); ++i) {
610
611 STATUS("processing: " << setw(3) << i->first << endl);
612
613 buffer_type data;
614
615 for (map<int, buffer_type>::iterator receiver = i->second.begin(); receiver != i->second.end(); ++receiver) {
616
617 // filter similar hits
618
619 sort(receiver->second.begin(), receiver->second.end(), JTransmission::compare(precision));
620
621 buffer_type::iterator __end = unique(receiver->second.begin(), receiver->second.end(), JTransmission::equals(precision));
622
623 // selection based on quality
624
625 for (buffer_type::const_iterator p = receiver->second.begin(); p != __end; ++p) {
626 if (p->getQ() >= parameters.Q * (parameters.Q <= 1.0 ? p->getW() : 1.0)) {
627 data.push_back(*p);
628 }
629 }
630 }
631
632 if (!data.empty()) {
633
634 sort(data.begin(), data.end(), make_comparator(&hit_type::getToA, JComparison::lt())); // sort according time-of-arrival
635
636 buffer_type buffer(factoryLimit); // local buffer of hits
637 event_type out[2]; // FIFO of events
638
639 TH1D* h1 = M1[i->first];
640 TH1D* h2 = M2[i->first];
641 TH1D* h3 = M3[i->first];
642 TH1D* h4 = M4[i->first];
643 TH1D* h5 = M5[i->first];
644
645 for (buffer_type::const_iterator p = data.begin(); p != data.end(); ++p) {
646
647 buffer_type::const_iterator q = p;
648
649 while (++q != data.end() && q->getToA() - p->getToA() <= parameters.TMax_s) {}
650
651 h1->Fill(distance(p,q));
652
653 if (distance(p,q) >= parameters.numberOfHits) {
654
655 if (distance(p,q) < factoryLimit) {
656
657 if ((int) out[1].size() >= parameters.numberOfHits && velo(out[1].getVertex(), p, q) >= parameters.numberOfHits) {
658
659 // ?
660
661 } else {
662
663 buffer_type::iterator root = buffer.begin();
664 buffer_type::iterator __p = buffer.begin();
665 buffer_type::iterator __q = buffer.begin();
666
667 *root = *p;
668
669 ++__p;
670 ++__q;
671
672 for (buffer_type::const_iterator i = p; ++i != q; ) {
673 if (match(*p,*i)) {
674 *__q = *i;
675 ++__q;
676 }
677 }
678
679 h2->Fill(distance(root,__q));
680
681 if (distance(root,__q) >= parameters.numberOfHits) {
682
683 __q = clusterize(__p, __q, match);
684
685 h3->Fill(distance(root,__q));
686
687 if (distance(root,__q) >= parameters.numberOfHits) {
688
689 vertex_type vx = Velo(center, *root, __p, __q, parameters.numberOfHits);
690
691 if (vx.N >= parameters.numberOfHits) {
692
693 vx = velo(vx.getPosition(), *root, __p, __q, parameters.numberOfHits);
694
695 h4->Fill(vx.N);
696
697 if (vx.N >= parameters.numberOfHits) {
698 out[1] = event_type(vx.getVertex(), JEvent(detector.getID(), out[0].getCounter() + 1, i->first, root, __q));
699 }
700 }
701 }
702 }
703 }
704
705 } else {
706
707 out[1] = event_type(JEvent(detector.getID(), out[0].getCounter() + 1, i->first, p, q));
708 }
709
710 if (!out[1].empty()) {
711
712 if (out[0].empty()) {
713
714 out[0] = out[1]; // shift
715
716 } else if (overlap2D(out[0],out[1])) {
717
718 out[0].merge(out[1]); // merge
719
720 } else {
721
722 STATUS("trigger: " << out[0] << endl);
723
724 h5->Fill(out[0].size());
725
726 queue.push_back(out[0]); // store
727
728 out[0] = out[1]; // shift
729 }
730
731 out[1].clear();
732 }
733 }
734 }
735
736 if (!out[0].empty()) {
737
738 queue.push_back(out[0]); // store
739 }
740 STATUS(endl);
741 }
742 }
743
744 if (!queue.empty()) {
745
746 sort(queue.begin(), queue.end());
747
748 for (vector<event_type>::iterator p = queue.begin(); p != queue.end(); ) {
749
750 // select event with highest quality within given time window
751
752 vector<event_type>::iterator q = p;
753 vector<event_type>::iterator i = p;
754
755 STATUS(endl);
756
757 for ( ; q != queue.end() && overlap1D(*p,*q); ++q) {
758
759 STATUS("| " << *q << endl);
760
761 if (q->getQ() > i->getQ()) {
762 i = q;
763 }
764 }
765
766 STATUS("event: " << *i << endl);
767
768 outputFile.put(*i);
769
770 G2[i->getID()].put(i->getX(), i->getY());
771
772 p = q;
773 }
774 }
775
776 JMultipleFileScanner<JMeta> io(inputFile);
777
778 io >> outputFile;
779
780 for (const auto M : { &M1, &M2, &M3, &M4, &M5 }){
781 for (const auto& i : *M) {
782 outputFile.put(*i.second);
783 }
784 }
785
786 for (map<int, JGraph_t>::const_iterator i = G2.begin(); i != G2.end(); ++i) {
787 outputFile.put(JGraph(i->second, MAKE_CSTRING("G[" << FILL(2,'0') << i->first << "]")));
788 }
789
790 outputFile.close();
791}
Acoustics support kit.
Acoustics toolkit.
int main(int argc, char **argv)
Constant sound velocity.
Acoustic event.
Acoustic receiver.
Sound velocity.
ROOT TTree parameter settings.
Acoustic event.
Acoustic transmission.
Acoustic trigger parameters.
3D match criterion for acoustic signals.
const JAbstractSoundVelocity & V
virtual bool operator()(const hit_type &first, const hit_type &second) const override
Match operator.
JMatch3D(const JAbstractSoundVelocity &V, const double Tmax_s=0.0)
Constructor.
Auxiliary classes and methods for acoustic position calibration.
JContainer< std::vector< JHydrophone > > hydrophones_container
Definition JSydney.cc:80
Interface for depth dependend velocity of sound.
virtual double getTime(const double D_m, const double z1, const double z2) const =0
Get propagation time of sound.
int getCounter() const
Get counter.
Definition JCounter.hh:50
Match of two events considering overlap in time and position.
JEventOverlap(const double Tmax_s, const double Dmax_m=std::numeric_limits< double >::max())
Constructor.
bool operator()(const event_type &first, const event_type &second) const
Match criterion.
Acoustic event.
Definition JEvent.hh:43
void merge(const JEvent &event)
Merge event.
Definition JEvent.hh:165
Acoustic receiver.
Definition JReceiver.hh:30
double getT(const double t_s) const
Get corrected time.
Definition JReceiver.hh:72
Implementation for depth dependend velocity of sound.
JSoundVelocity & set(const double z0)
Set depth.
Time-of-arrival data from acoustic piezo sensor or hydrophone.
Definition JToA.hh:28
uint32_t DOMID
DAQ run number.
Definition JToA.hh:34
uint32_t QUALITYFACTOR
The ticks (16ns) part of the DAQ frame timestamp.
Definition JToA.hh:39
uint32_t QUALITYNORMALISATION
A measure of how good the waveform match was to the signal.
Definition JToA.hh:40
int32_t WAVEFORMID
DOM unique identifeir.
Definition JToA.hh:35
int32_t DETID
Definition JToA.hh:32
int32_t RUN
detector identifier
Definition JToA.hh:33
double TOA_S() const
Time of Arrival, expressed in seconds relative to Unix epoch (1 January 1970 00:00:00 UTC)
Definition JToAImp.cc:40
Auxiliary class to compare transmissions.
Auxiliary class to compare transmissions.
Acoustic transmission.
double getToA() const
Get calibrated time of arrival.
double getQ() const
Get quality.
double Q
minimal quality if larger than one; else minimal normalised quality
double TMax_s
maximal difference between times of emission [s]
int numberOfHits
minimal number of hits to trigger event
int operator()(const JVertex3D &vx, T __begin, T __end) const
Check vertex.
vertex_type operator()(const JPosition3D &position, const hit_type &root, T __begin, T __end, const int numberOfHits=0) const
Locate vertex around given position.
const JAbstractSoundVelocity & V
JVelo(const JAbstractSoundVelocity &V, const double RMax_m, const double Xv_m, const double factor=1.0)
Constructor.
double getQ() const
Get average quality.
event_type(const JEvent &event)
Copy constructor.
friend std::ostream & operator<<(std::ostream &out, const event_type &object)
Write event to output stream.
static double TOA_s
start time of data
event_type(const JVertex3D &vertex, const JEvent &event)
Constructor.
Transmission with position.
hit_type(const JPosition3D &position, const JTransmission &transmission)
Constructor.
vertex_type(const JVector3D &p0, const double t0, const int n0, const double q0)
Constructor.