Jpp  18.2.0-rc.1
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
JHydrophone.cc File Reference

Example program to plot hydrophone data. More...

#include <iostream>
#include <iomanip>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "JDB/JToAshort.hh"
#include "JDB/JSupport.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JTripod.hh"
#include "JDetector/JTransmitter.hh"
#include "JDetector/JHydrophone.hh"
#include "JSupport/JMultipleFileScanner.hh"
#include "JSupport/JTreeScanner.hh"
#include "JTools/JHashMap.hh"
#include "JTools/JRange.hh"
#include "JLang/JComparator.hh"
#include "JLang/JComparison.hh"
#include "JROOT/JManager.hh"
#include "JROOT/JRootToolkit.hh"
#include "JAcoustics/JEmitterID.hh"
#include "JAcoustics/JSoundVelocity.hh"
#include "JAcoustics/JEmitter.hh"
#include "JAcoustics/JReceiver.hh"
#include "JAcoustics/JTransceiver.hh"
#include "JAcoustics/JAcousticsToolkit.hh"
#include "JAcoustics/JAcousticsSupportkit.hh"
#include "JAcoustics/JHit.hh"
#include "JAcoustics/JEvent.hh"
#include "JAcoustics/JSupport.hh"
#include "Jeep/JContainer.hh"
#include "Jeep/JPrint.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 plot hydrophone data.

Author
mdejong

Definition in file JHydrophone.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 120 of file JHydrophone.cc.

121 {
122  using namespace std;
123  using namespace JPP;
124 
128 
130  string detectorFile;
131  JLimit_t& numberOfEvents = inputFile.getLimit();
132  string outputFile;
133  JSoundVelocity V = getSoundVelocity; // default sound velocity
134  tripods_container tripods; // tripods
135  transmitters_container transmitters; // transmitters
136  hydrophones_container hydrophones; // hydrophones
137  double Q;
138  int debug;
139 
140  try {
141 
142  JParser<> zap("Example program to plot hydrophone data.");
143 
144  zap['f'] = make_field(inputFile);
145  zap['n'] = make_field(numberOfEvents) = JLimit::max();
146  zap['o'] = make_field(outputFile) = "hydrophone.root";
147  zap['a'] = make_field(detectorFile);
148  zap['V'] = make_field(V, "sound velocity") = JPARSER::initialised();
149  zap['T'] = make_field(tripods, "tripod data");
150  zap['Y'] = make_field(transmitters, "transmitter data") = JPARSER::initialised();
151  zap['H'] = make_field(hydrophones, "hydrophone data") = JPARSER::initialised();
152  zap['W'] = make_field(getEmitterID, "waveform identification data") = JPARSER::initialised();
153  zap['Q'] = make_field(Q, "quality") = 0.0;
154  zap['d'] = make_field(debug) = 2;
155 
156  zap(argc, argv);
157  }
158  catch(const exception &error) {
159  FATAL(error.what() << endl);
160  }
161 
162 
164 
165  try {
166  load(detectorFile, detector);
167  }
168  catch(const JException& error) {
169  FATAL(error);
170  }
171 
172  V.set(detector.getUTMZ());
173 
174  JHashMap<int, JReceiver> receivers;
175  JHashMap<int, JEmitter> emitters;
176 
177  for (hydrophones_container::const_iterator i = hydrophones.begin(); i != hydrophones.end(); ++i) {
178 
179  for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) {
180 
181  if (i->getLocation() == module->getLocation()) {
182 
183  receivers[module->getID()] = JReceiver(module->getID(),
184  module->getPosition() + i->getPosition(),
185  module->getT0() * 1.0e-9);
186  break;
187  }
188  }
189  }
190 
191  for (tripods_container::const_iterator i = tripods.begin(); i != tripods.end(); ++i) {
192  emitters[i->getID()] = JEmitter(i->getID(),
193  i->getUTMPosition() - detector.getUTMPosition());
194  }
195 
196  for (transmitters_container::const_iterator i = transmitters.begin(); i != transmitters.end(); ++i) {
197  try {
198  emitters[i->getID()] = JEmitter(i->getID(),
199  i->getPosition() + detector.getModule(i->getLocation()).getPosition());
200  }
201  catch(const exception&) {
202  continue; // if no string available, discard transmitter
203  }
204  }
205 
206 
207  const JRange<double> T0(-0.05, +0.05);
208  const JRange<double> T1(-0.01, +0.01);
209 
210  TH1D h0(MAKE_CSTRING("Q0 " << T0), NULL, 100, 0.0, 8.0);
211  TH1D h1(MAKE_CSTRING("Q1 " << T1), NULL, 100, 0.0, 8.0);
212 
213  JManager<key_type, TH1D> H1(new TH1D("[%]", NULL, 50000, -5.0e-2, +5.0e-2));
214 
215 
216  while (inputFile.hasNext()) {
217 
218  STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl);
219 
220  const JEvent* evt = inputFile.next();
221 
222  if (!evt->empty() && emitters.has(evt->getID())) {
223 
225 
226  for (JEvent::const_iterator hit = evt->begin(); hit != evt->end(); ++hit) {
227  if (receivers.has(hit->getID())) {
228  buffer[hit->getID()].push_back(*hit);
229  }
230  }
231 
232  for (map<int, vector<JTransmission> >::iterator i = buffer.begin(); i != buffer.end(); ++i) {
233 
234  sort(i->second.begin(), i->second.end(), make_comparator(&JTransmission::getQ, JComparison::gt()));
235 
236  vector<JTransmission>::const_iterator hit = i->second.begin();
237 
238  if (hit->getQ() >= Q) {
239 
240  const JPosition3D& p1 = emitters [evt->getID()].getPosition();
241  const JPosition3D& p2 = receivers[hit->getID()].getPosition();
242 
243  const double D = p2.getDistance(p1);
244  const double Vi = V.getInverseVelocity(D, p1.getZ(), p2.getZ());
245 
246  const double t0 = evt->begin()->getToE() + D * Vi;
247  const double t1 = hit->getToA() - t0;
248 
249  H1[key_type(hit->getID(), evt->getID())]->Fill(t1);
250 
251  const double Q = log10(hit->getQ());
252 
253  if (T0(t1)) { h0.Fill(Q); }
254  if (T1(t1)) { h1.Fill(Q); }
255  }
256  }
257  }
258  }
259  STATUS(endl);
260 
261  TFile out(outputFile.c_str(), "recreate");
262 
263  out << h0 << h1 << H1;
264 
265  out.Write();
266  out.Close();
267 }
Utility class to parse command line options.
Definition: JParser.hh:1514
General exception.
Definition: JException.hh:24
Acoustic receiver.
Definition: JReceiver.hh:27
Q(UTCMax_s-UTCMin_s)-livetime_s
TPaveText * p1
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.
JContainer< std::vector< JTransmitter > > transmitters_container
Definition: JSydney.cc:79
General purpose class for hash map of unique keys.
Definition: JHashMap.hh:72
#define STATUS(A)
Definition: JMessage.hh:63
Detector data structure.
Definition: JDetector.hh:89
#define MAKE_CSTRING(A)
Make C-string.
Definition: JPrint.hh:136
then fatal Number of tripods
Definition: JFootprint.sh:45
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition: JParser.hh:83
V(JDAQEvent-JTriggerReprocessor)*1.0/(JDAQEvent+1.0e-10)
string outputFile
double getDistance(const JVector3D &pos) const
Get distance to point.
Definition: JVector3D.hh:270
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.
Auxiliary class to manage set of compatible ROOT objects (e.g. histograms) using unique keys...
Definition: JManager.hh:43
Detector file.
Definition: JHead.hh:226
JContainer< std::vector< JHydrophone > > hydrophones_container
Definition: JSydney.cc:78
Acoustic emitter.
Definition: JEmitter.hh:27
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:1989
set_variable E_E log10(E_{fit}/E_{#mu})"
JPosition3D getPosition(const Vec &pos)
Get position.
JContainer< std::vector< JTripod > > tripods_container
Definition: JSydney.cc:77
double getQ(const double D_m, const double f_kHz, const double d_m)
Get relative quality for given frequency at given distance.
Implementation for depth dependend velocity of sound.
#define FATAL(A)
Definition: JMessage.hh:67
p2
Definition: module-Z:fit.sh:74
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
General purpose class for object reading from a list of file names.
Acoustic event.
Data structure for position in three dimensions.
Definition: JPosition3D.hh:36
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:84
static JEmitterID getEmitterID
Function object for emitter identification.
Definition: JEmitterID.hh:119
do set_variable DETECTOR_TXT $WORKDIR detector
int getID() const
Get emitter identifier.
do echo Generating $dir eval D
Definition: JDrawLED.sh:53
double getZ() const
Get z position.
Definition: JVector3D.hh:115
int debug
debug level
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62