Jpp  19.0.0
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 "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/JToA.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 118 of file JHydrophone.cc.

119 {
120  using namespace std;
121  using namespace JPP;
122 
126 
128  string detectorFile;
129  JLimit_t& numberOfEvents = inputFile.getLimit();
130  string outputFile;
131  JSoundVelocity V = getSoundVelocity; // default sound velocity
132  tripods_container tripods; // tripods
133  transmitters_container transmitters; // transmitters
134  hydrophones_container hydrophones; // hydrophones
135  double Q;
136  int debug;
137 
138  try {
139 
140  JParser<> zap("Example program to plot hydrophone data.");
141 
142  zap['f'] = make_field(inputFile);
143  zap['n'] = make_field(numberOfEvents) = JLimit::max();
144  zap['o'] = make_field(outputFile) = "hydrophone.root";
145  zap['a'] = make_field(detectorFile);
146  zap['V'] = make_field(V, "sound velocity") = JPARSER::initialised();
147  zap['T'] = make_field(tripods, "tripod data");
148  zap['Y'] = make_field(transmitters, "transmitter data") = JPARSER::initialised();
149  zap['H'] = make_field(hydrophones, "hydrophone data") = JPARSER::initialised();
150  zap['W'] = make_field(getEmitterID, "waveform identification data") = JPARSER::initialised();
151  zap['Q'] = make_field(Q, "quality") = 0.0;
152  zap['d'] = make_field(debug) = 2;
153 
154  zap(argc, argv);
155  }
156  catch(const exception &error) {
157  FATAL(error.what() << endl);
158  }
159 
160 
162 
163  try {
164  load(detectorFile, detector);
165  }
166  catch(const JException& error) {
167  FATAL(error);
168  }
169 
170  V.set(detector.getUTMZ());
171 
172  JHashMap<int, JReceiver> receivers;
173  JHashMap<int, JEmitter> emitters;
174 
175  for (hydrophones_container::const_iterator i = hydrophones.begin(); i != hydrophones.end(); ++i) {
176 
177  for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) {
178 
179  if (i->getLocation() == module->getLocation()) {
180 
181  receivers[module->getID()] = JReceiver(module->getID(),
182  module->getPosition() + i->getPosition(),
183  module->getT0() * 1.0e-9);
184  break;
185  }
186  }
187  }
188 
189  for (tripods_container::const_iterator i = tripods.begin(); i != tripods.end(); ++i) {
190  emitters[i->getID()] = JEmitter(i->getID(),
191  i->getUTMPosition() - detector.getUTMPosition());
192  }
193 
194  for (transmitters_container::const_iterator i = transmitters.begin(); i != transmitters.end(); ++i) {
195  try {
196  emitters[i->getID()] = JEmitter(i->getID(),
197  i->getPosition() + detector.getModule(i->getLocation()).getPosition());
198  }
199  catch(const exception&) {
200  continue; // if no string available, discard transmitter
201  }
202  }
203 
204 
205  const JRange<double> T0(-0.05, +0.05);
206  const JRange<double> T1(-0.01, +0.01);
207 
208  TH1D h0(MAKE_CSTRING("Q0 " << T0), NULL, 100, 0.0, 8.0);
209  TH1D h1(MAKE_CSTRING("Q1 " << T1), NULL, 100, 0.0, 8.0);
210 
211  JManager<key_type, TH1D> H1(new TH1D("[%]", NULL, 50000, -5.0e-2, +5.0e-2));
212 
213 
214  while (inputFile.hasNext()) {
215 
216  STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl);
217 
218  const JEvent* evt = inputFile.next();
219 
220  if (!evt->empty() && emitters.has(evt->getID())) {
221 
223 
224  for (JEvent::const_iterator hit = evt->begin(); hit != evt->end(); ++hit) {
225  if (receivers.has(hit->getID())) {
226  buffer[hit->getID()].push_back(*hit);
227  }
228  }
229 
230  for (map<int, vector<JTransmission> >::iterator i = buffer.begin(); i != buffer.end(); ++i) {
231 
232  sort(i->second.begin(), i->second.end(), make_comparator(&JTransmission::getQ, JComparison::gt()));
233 
234  vector<JTransmission>::const_iterator hit = i->second.begin();
235 
236  if (hit->getQ() >= Q) {
237 
238  const JPosition3D& p1 = emitters [evt->getID()].getPosition();
239  const JPosition3D& p2 = receivers[hit->getID()].getPosition();
240 
241  const double D = p2.getDistance(p1);
242  const double Vi = V.getInverseVelocity(D, p1.getZ(), p2.getZ());
243 
244  const double t0 = evt->begin()->getToE() + D * Vi;
245  const double t1 = hit->getToA() - t0;
246 
247  H1[key_type(hit->getID(), evt->getID())]->Fill(t1);
248 
249  const double Q = log10(hit->getQ());
250 
251  if (T0(t1)) { h0.Fill(Q); }
252  if (T1(t1)) { h1.Fill(Q); }
253  }
254  }
255  }
256  }
257  STATUS(endl);
258 
259  TFile out(outputFile.c_str(), "recreate");
260 
261  out << h0 << h1 << H1;
262 
263  out.Write();
264  out.Close();
265 }
Utility class to parse command line options.
Definition: JParser.hh:1711
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:80
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:84
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:79
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:2158
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:78
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