Jpp  17.2.0
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
JFremantle.cc File Reference

Application to make a global fit of the detector geometry to acoustic data. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <vector>
#include <set>
#include <map>
#include <deque>
#include <algorithm>
#include <type_traits>
#include <functional>
#include <future>
#include <mutex>
#include <thread>
#include <queue>
#include "TROOT.h"
#include "TFile.h"
#include "JLang/JPredicate.hh"
#include "JLang/JComparator.hh"
#include "JLang/JComparison.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JTripod.hh"
#include "JDetector/JTransmitter.hh"
#include "JDetector/JModule.hh"
#include "JDetector/JHydrophone.hh"
#include "JTools/JHashMap.hh"
#include "JTools/JRange.hh"
#include "JSupport/JMultipleFileScanner.hh"
#include "JSupport/JFileRecorder.hh"
#include "JSupport/JMeta.hh"
#include "JAcoustics/JSoundVelocity.hh"
#include "JAcoustics/JEmitter.hh"
#include "JAcoustics/JAcousticsToolkit.hh"
#include "JAcoustics/JHit.hh"
#include "JAcoustics/JFitParameters.hh"
#include "JAcoustics/JKatoomba.hh"
#include "JAcoustics/JEvent.hh"
#include "JAcoustics/JEvt.hh"
#include "JAcoustics/JEvtToolkit.hh"
#include "JAcoustics/JSupport.hh"
#include "JAcoustics/JTransmission_t.hh"
#include "Jeep/JContainer.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

Application to make a global fit of the detector geometry to acoustic data.


Author
mdejong

Definition in file JFremantle.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 221 of file JFremantle.cc.

222 {
223  using namespace std;
224  using namespace JPP;
225 
226  typedef JContainer< vector<JTripod> > tripods_container;
227  typedef JContainer< vector<JTransmitter> > transmitters_container;
228  typedef JContainer< vector<JHydrophone> > hydrophones_container;
229  typedef JContainer< set<JTransmission_t> > disable_container;
230 
233  string detectorFile;
234  JLimit_t& numberOfEvents = inputFile.getLimit();
235  JSoundVelocity V = getSoundVelocity; // default sound velocity
236  tripods_container tripods; // tripods
237  transmitters_container transmitters; // transmitters
238  hydrophones_container hydrophones; // hydrophones
239  JFitParameters parameters; // fit parameters
240  bool unify; // unify weighing of pings
241  disable_container disable; // disable tansmissions
242  size_t jobs; // number of parallel jobs
243  int sleep_us; // sleep time [us]
244  int debug;
245 
246  try {
247 
248  JParser<> zap("Application to fit position calibration model to acoustic data.");
249 
250  zap['f'] = make_field(inputFile, "output of JAcousticEventBuilder[.sh]");
251  zap['o'] = make_field(outputFile);
252  zap['n'] = make_field(numberOfEvents) = JLimit::max();
253  zap['a'] = make_field(detectorFile);
255  zap['V'] = make_field(V, "sound velocity") = JPARSER::initialised();
256  zap['T'] = make_field(tripods, "tripod data");
257  zap['Y'] = make_field(transmitters, "transmitter data") = JPARSER::initialised();
258  zap['H'] = make_field(hydrophones, "hydrophone data") = JPARSER::initialised();
259  zap['M'] = make_field(getMechanics, "mechanics data") = JPARSER::initialised();
260  zap['u'] = make_field(unify, "unify weighing of pings");
261  zap['!'] = make_field(disable, "disable transmission") = JPARSER::initialised();
262  zap['N'] = make_field(jobs, "number of parallel jobs") = 1;
263  zap['s'] = make_field(sleep_us, "sleep time [us]") = 100;
264  zap['d'] = make_field(debug) = 1;
265 
266  zap(argc, argv);
267  }
268  catch(const exception &error) {
269  FATAL(error.what() << endl);
270  }
271 
272  ROOT::EnableThreadSafety();
273 
274  cout.tie(&cerr);
275 
277 
278  try {
279  load(detectorFile, detector);
280  }
281  catch(const JException& error) {
282  FATAL(error);
283  }
284 
285  JHashMap<int, JLocation> receivers;
286  JHashMap<int, JEmitter> emitters;
287 
288  for (JDetector::const_iterator i = detector.begin(); i != detector.end(); ++i) {
289  receivers[i->getID()] = i->getLocation();
290  }
291 
292  for (tripods_container::const_iterator i = tripods.begin(); i != tripods.end(); ++i) {
293  emitters[i->getID()] = JEmitter(i->getID(),
294  i->getUTMPosition() - detector.getUTMPosition());
295  }
296 
297  for (transmitters_container::const_iterator i = transmitters.begin(); i != transmitters.end(); ++i) {
298  try {
299  emitters[i->getID()] = JEmitter(i->getID(),
300  i->getPosition() + detector.getModule(i->getLocation()).getPosition());
301  }
302  catch(const exception&) {
303  continue; // if no module available, discard transmitter
304  }
305  }
306 
307  V.set(detector.getUTMZ()); // sound velocity at detector depth
308 
309  JGeometry geometry(detector, hydrophones);
310  JKatoomba_t katoomba(geometry, V, parameters);
311 
314 
317 
318  DEBUG(geometry);
319 
320  string oid; // detector identifier
321 
322  { // sort input files
323  map<double, string> zmap;
324 
325  for (const string& file_name : inputFile) {
326 
327  STATUS(file_name << '\r'); DEBUG(endl);
328 
329  for (JMultipleFileScanner<JEvent> in(file_name, 1); in.hasNext(); ) {
330 
331  const JEvent* evt = in.next();
332 
333  if (oid == "")
334  oid = evt->getOID();
335  else if (oid != evt->getOID()) // consistency check
336  FATAL("Invalid detector identifier " << evt->getOID() << " != " << oid << endl);
337 
338  if (!evt->empty()) {
339  zmap[evt->begin()->getToE()] = file_name;
340  }
341  }
342  }
343  STATUS(endl);
344 
345  inputFile.clear();
346 
347  for (map<double, string>::const_iterator i = zmap.begin(); i != zmap.end(); ++i) {
348  inputFile.push_back(i->second);
349  }
350  }
351 
352  outputFile.open();
353 
354  outputFile.put(JMeta(argc, argv));
355  outputFile.put(parameters);
356 
357  try {
358  JFremantle fremantle(oid, katoomba, outputFile, jobs);
359 
360  typedef deque<JEvent> buffer_type;
361 
362  for (buffer_type zbuf; inputFile.hasNext(); ) {
363 
364  STATUS(inputFile.getFilename() << '\r'); DEBUG(endl);
365 
366  // read one file at a time
367 
368  for (const string file_name = inputFile.getFilename(); inputFile.hasNext() && file_name == inputFile.getFilename(); ) {
369 
370  const JEvent* evt = inputFile.next();
371 
372  if (emitters.has(evt->getID())) {
373  zbuf.push_back(*evt);
374  }
375  }
376 
377  sort(zbuf.begin(), zbuf.end()); // sort according first time-of-emission
378 
379  for (buffer_type::iterator p = zbuf.begin(), q; p != zbuf.end(); p = q) {
380 
381  for (q = p; ++q != zbuf.end() && q->begin()->getToE() <= p->rbegin()->getToE() + parameters.Tmax_s; ) {}
382 
383  if (q == zbuf.end()) {
384 
385  if (inputFile.hasNext()) {
386 
387  zbuf.erase(zbuf.begin(), p); // remove processed data and continue reading
388 
389  break;
390  }
391  }
392 
393  if (getNumberOfEmitters(p,q) >= parameters.Nmin) {
394 
395  map<int, int> numberOfPings;
396 
397  for (buffer_type::const_iterator i = p; i != q; ++i) {
398  numberOfPings[i->getID()] += 1;
399  }
400 
401  for (map<int, int>::const_iterator i = numberOfPings.begin(); i != numberOfPings.end(); ++i) {
402  DEBUG("Number of pings " << setw(2) << i->first << ' ' << setw(3) << i->second << endl);
403  }
404 
405  int minimum_number_of_pings = numeric_limits<int>::max();
406 
407  for (map<int, int>::const_iterator i = numberOfPings.begin(); i != numberOfPings.end(); ++i) {
408  minimum_number_of_pings = min(minimum_number_of_pings, i->second);
409  }
410 
411  data_type data;
412 
413  for (buffer_type::iterator evt = p; evt != q; ++evt) {
414 
415  sort(evt->begin(), evt->end(), compare);
416 
417  JEvent::iterator __end = unique(evt->begin(), evt->end(), make_comparator(&JTransmission::getID, JComparison::eq()));
418 
419  const JEmitter& emitter = emitters[evt->getID()];
420  const double signal = (unify ? (double) minimum_number_of_pings / (double) numberOfPings[evt->getID()] : 1.0);
421 
422  for (JEvent::const_iterator i = evt->begin(); i != __end; ++i) {
423 
424  if (disable.count(JTransmission_t(evt->getID(), i->getID())) == 0 &&
425  disable.count(JTransmission_t(-1, i->getID())) == 0) {
426 
427  if (receivers.has(i->getID()) && geometry.hasLocation(receivers[i->getID()]) && i->getQ() >= parameters.Qmin) {
428 
429  data.push_back(hit_type(emitter,
430  distance(zbuf.begin(),evt),
431  receivers[i->getID()],
432  JPDFGauss(i->getToA(), parameters.sigma_s, signal, parameters.background)));
433  }
434  }
435  }
436  }
437 
438  while (fremantle.backlog() > jobs) {
439  this_thread::sleep_for(chrono::microseconds(sleep_us));
440  }
441 
442  fremantle.enqueue(data);
443  }
444  }
445  }
446  STATUS(endl);
447  }
448  catch(const exception& error) {
449  FATAL("main " << error.what());
450  }
451 
452  JMultipleFileScanner<JMeta> io(inputFile);
453 
454  io >> outputFile;
455 
456  outputFile.close();
457 }
Worker class for fit function of acoustic model.
Definition: JKatoomba.hh:869
Auxiliary class for ROOT I/O of application specific meta data.
Definition: JMeta.hh:70
Object writing to file.
Utility class to parse command line options.
Definition: JParser.hh:1517
General exception.
Definition: JException.hh:23
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
static JDetectorMechanics getMechanics
Function object to get string mechanics.
Definition: JMechanics.hh:243
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
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
*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
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
size_t getNumberOfEmitters(T __begin, T __end)
Get number of emitters.
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.
Detector file.
Definition: JHead.hh:226
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:1993
const std::string & getOID() const
Get detector identifier.
JPosition3D getPosition(const Vec &pos)
Get position.
static struct JACOUSTICS::@4 compare
Auxiliary data structure to sort transmissions.
Implementation for depth dependend velocity of sound.
#define FATAL(A)
Definition: JMessage.hh:67
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.
Custom probability density function of time-of-arrival.
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:73
do echo n Creating graphics for string $STRING for((FLOOR=$FIRST_FLOOR;$FLOOR<=$LAST_FLOOR;FLOOR+=1))
do set_variable DETECTOR_TXT $WORKDIR detector
int getID() const
Get identifier.
then fatal Wrong number of arguments fi set_variable DETECTOR $argv[1] set_variable INPUT_FILE $argv[2] eval JPrintDetector a $DETECTOR O IDENTIFIER eval JPrintDetector a $DETECTOR O SUMMARY JAcoustics sh $DETECTOR_ID source JAcousticsToolkit sh CHECK_EXIT_CODE typeset A EMITTERS get_tripods $WORKDIR tripod txt EMITTERS get_transmitters $WORKDIR transmitter txt EMITTERS for EMITTER in
Definition: JCanberra.sh:46
Acoustic transmission identifier.
Template definition of fit function of acoustic model.
Definition: JKatoomba.hh:115
int debug
debug level
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62