Jpp  17.1.1
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
software/JAcoustics/JKatoomba.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 "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TH2D.h"
#include "JLang/JPredicate.hh"
#include "JLang/JComparator.hh"
#include "JLang/JComparison.hh"
#include "JROOT/JManager.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 "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 software/JAcoustics/JKatoomba.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 79 of file software/JAcoustics/JKatoomba.cc.

80 {
81  using namespace std;
82  using namespace JPP;
83 
84  typedef JContainer< vector<JTripod> > tripods_container;
85  typedef JContainer< vector<JTransmitter> > transmitters_container;
86  typedef JContainer< vector<JHydrophone> > hydrophones_container;
87  typedef JContainer< set<JTransmission_t> > disable_container;
88 
91  string detectorFile;
92  JLimit_t& numberOfEvents = inputFile.getLimit();
93  JSoundVelocity V = getSoundVelocity; // default sound velocity
94  tripods_container tripods; // tripods
95  transmitters_container transmitters; // transmitters
96  hydrophones_container hydrophones; // hydrophones
97  JFitParameters parameters; // fit parameters
98  int fit; // type of fit
99  bool unify; // unify weighing of pings
100  disable_container disable; // disable tansmissions
101  int debug;
102 
103  try {
104 
105  JParser<> zap("Application to fit position calibration model to acoustic data.");
106 
107  zap['f'] = make_field(inputFile, "output of JAcousticEventBuilder[.sh]");
108  zap['o'] = make_field(outputFile);
109  zap['n'] = make_field(numberOfEvents) = JLimit::max();
110  zap['a'] = make_field(detectorFile);
112  zap['V'] = make_field(V, "sound velocity") = JPARSER::initialised();
113  zap['T'] = make_field(tripods, "tripod data");
114  zap['Y'] = make_field(transmitters, "transmitter data") = JPARSER::initialised();
115  zap['H'] = make_field(hydrophones, "hydrophone data") = JPARSER::initialised();
116  zap['M'] = make_field(getMechanics, "mechanics data") = JPARSER::initialised();
117  zap['F'] = make_field(fit, "fit type") = linear_t, simplex_t, gandalf_t;
118  zap['u'] = make_field(unify, "unify weighing of pings");
119  zap['!'] = make_field(disable, "disable transmission") = JPARSER::initialised();
120  zap['d'] = make_field(debug) = 1;
121 
122  zap(argc, argv);
123  }
124  catch(const exception &error) {
125  FATAL(error.what() << endl);
126  }
127 
128 
129  cout.tie(&cerr);
130 
132 
133  try {
134  load(detectorFile, detector);
135  }
136  catch(const JException& error) {
137  FATAL(error);
138  }
139 
141 
142  JHashMap<int, JLocation> receivers;
143  JHashMap<int, JEmitter> emitters;
144 
145  for (JDetector::const_iterator i = detector.begin(); i != detector.end(); ++i) {
146  receivers[i->getID()] = i->getLocation();
147  }
148 
149  for (tripods_container::const_iterator i = tripods.begin(); i != tripods.end(); ++i) {
150  emitters[i->getID()] = JEmitter(i->getID(),
151  i->getUTMPosition() - detector.getUTMPosition());
152  }
153 
154  for (transmitters_container::const_iterator i = transmitters.begin(); i != transmitters.end(); ++i) {
155  try {
156  emitters[i->getID()] = JEmitter(i->getID(),
157  i->getPosition() + detector.getModule(i->getLocation()).getPosition());
158  }
159  catch(const exception&) {
160  continue; // if no module available, discard transmitter
161  }
162  }
163 
164  V.set(detector.getUTMZ()); // sound velocity at detector depth
165 
166  JGeometry geometry(detector, hydrophones);
167 
168  DEBUG(geometry);
169 
170  JKatoomba<JEstimator> estimator(geometry, V, parameters.option);
171  JKatoomba<JAbstractMinimiser> evaluator(geometry, V, parameters.option);
172  JKatoomba<JSimplex> simplex (geometry, V, parameters.option);
173  JKatoomba<JGandalf> gandalf (geometry, V, parameters.option);
174 
175  evaluator.estimator.reset(getMEstimator(parameters.mestimator));
176  simplex .estimator.reset(getMEstimator(parameters.mestimator));
177  gandalf .estimator.reset(getMEstimator(parameters.mestimator));
178 
179  simplex.debug = debug;
180  gandalf.debug = debug;
181 
182  JSimplex <JModel> ::MAXIMUM_ITERATIONS = 10000;
184 
185  typedef JHit<JPDFGauss> hit_type;
186 
187  TH1D h0("chi2/NDF", NULL, 50000, 0.0, 1000.0);
188  TH1D h1("h1", NULL, 51, -0.5, 50.5);
189  TH1D hn("hn", NULL, 100, 0.0, 6.0);
190 
191  JManager<int, TH2D> HA(new TH2D("ha[%]", NULL,
193  range.getLength() + 1, range.getLowerLimit() - 0.5, range.getUpperLimit() + 0.5));
194 
195  JManager<int, TH2D> HB(new TH2D("hb[%]", NULL,
197  range.getLength() + 1, range.getLowerLimit() - 0.5, range.getUpperLimit() + 0.5));
198 
199  for (Int_t i = 1; i <= HA->GetXaxis()->GetNbins(); ++i) {
200  HA->GetXaxis()->SetBinLabel(i, MAKE_CSTRING(geometry.at(i-1).first));
201  HB->GetXaxis()->SetBinLabel(i, MAKE_CSTRING(geometry.at(i-1).first));
202  }
203 
204  string oid; // detector identifier
205 
206  { // sort input files
207  map<double, string> zmap;
208 
209  for (const string& file_name : inputFile) {
210 
211  STATUS(file_name << '\r'); DEBUG(endl);
212 
213  for (JMultipleFileScanner<JEvent> in(file_name, 1); in.hasNext(); ) {
214 
215  const JEvent* evt = in.next();
216 
217  if (oid == "")
218  oid = evt->getOID();
219  else if (oid != evt->getOID()) // consistency check
220  FATAL("Invalid detector identifier " << evt->getOID() << " != " << oid << endl);
221 
222  if (!evt->empty()) {
223  zmap[evt->begin()->getToE()] = file_name;
224  }
225  }
226  }
227  STATUS(endl);
228 
229  inputFile.clear();
230 
231  for (map<double, string>::const_iterator i = zmap.begin(); i != zmap.end(); ++i) {
232  inputFile.push_back(i->second);
233  }
234  }
235 
236 
237  outputFile.open();
238 
239  outputFile.put(JMeta(argc, argv));
240  outputFile.put(parameters);
241 
242  int counter[] = { 0, 0 };
243 
244  typedef deque<JEvent> buffer_type;
245 
246  for (buffer_type zbuf; inputFile.hasNext(); ) {
247 
248  STATUS(inputFile.getFilename() << '\r'); DEBUG(endl);
249 
250  // read one file at a time
251 
252  for (const string file_name = inputFile.getFilename(); inputFile.hasNext() && file_name == inputFile.getFilename(); ) {
253 
254  const JEvent* evt = inputFile.next();
255 
256  if (emitters.has(evt->getID())) {
257  zbuf.push_back(*evt);
258  }
259  }
260 
261  sort(zbuf.begin(), zbuf.end()); // sort according first time-of-emission
262 
263  for (buffer_type::iterator p = zbuf.begin(), q; p != zbuf.end(); p = q) {
264 
265  for (q = p; ++q != zbuf.end() && q->begin()->getToE() <= p->rbegin()->getToE() + parameters.Tmax_s; ) {}
266 
267  if (q == zbuf.end()) {
268 
269  if (inputFile.hasNext()) {
270 
271  zbuf.erase(zbuf.begin(), p); // remove processed data and continue reading
272 
273  break;
274  }
275  }
276 
277  if (getNumberOfEmitters(p,q) >= parameters.Nmin) {
278 
279  h1.Fill((double) getNumberOfEmitters(p,q));
280 
281  map<int, int> numberOfPings;
282 
283  for (buffer_type::const_iterator i = p; i != q; ++i) {
284  numberOfPings[i->getID()] += 1;
285  }
286 
287  for (map<int, int>::const_iterator i = numberOfPings.begin(); i != numberOfPings.end(); ++i) {
288  DEBUG("Number of pings " << setw(2) << i->first << ' ' << setw(3) << i->second << endl);
289  }
290 
291  int minimum_number_of_pings = numeric_limits<int>::max();
292 
293  for (map<int, int>::const_iterator i = numberOfPings.begin(); i != numberOfPings.end(); ++i) {
294  minimum_number_of_pings = min(minimum_number_of_pings, i->second);
295  }
296 
297  vector<hit_type> data;
298 
299  for (buffer_type::iterator i = p; i != q; ++i) {
300 
301  {
302  const JEmitter& emitter = emitters[i->getID()];
303 
304  const double signal = (unify ? (double) minimum_number_of_pings / (double) numberOfPings[i->getID()] : 1.0);
305 
306  // select transmission with highest quality
307 
308  map<int, vector<JTransmission> > buffer;
309 
310  for (JEvent::const_iterator hit = i->begin(); hit != i->end(); ++hit) {
311  if (disable.count(JTransmission_t(i->getID(), hit->getID())) == 0 &&
312  disable.count(JTransmission_t(-1, hit->getID())) == 0) {
313  buffer[hit->getID()].push_back(*hit);
314  }
315  }
316 
317  for (map<int, vector<JTransmission> >::iterator ps = buffer.begin(); ps != buffer.end(); ++ps) {
318 
319  if (receivers.has(ps->first) && geometry.hasLocation(receivers[ps->first])) {
320 
321  sort(ps->second.begin(), ps->second.end(), make_comparator(&JTransmission::getQ, JComparison::gt()));
322 
323  if (ps->second.begin()->getQ() >= parameters.Qmin) {
324 
325  data.push_back(hit_type(emitter,
326  distance(p,i),
327  receivers[ps->first],
328  JPDFGauss(ps->second.begin()->getToA(), parameters.sigma_s, signal, parameters.background)));
329  }
330  }
331  }
332  }
333  }
334 
335 
336  JModel result = estimator(data.begin(), data.end()); // pre-fit
337 
338  double chi2 = evaluator(result, data.begin(), data.end());
339  double ndf = getWeight(data.begin(), data.end()) - result.getN();
340 
341  DEBUG("prefit:" << ' '
342  << FIXED(12,3) << chi2 << '/' << FIXED(6,1) << ndf << endl
343  << result << endl);
344 
345  for (vector<hit_type>::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
346  DEBUG("hit: " << *hit << ' ' << FIXED(9,3) << evaluator(result, *hit) << endl);
347  }
348 
349  for (vector<hit_type>::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
350  HA[hit->getID()]->Fill(geometry.getIndex(hit->getString()), hit->getFloor(), 1.0);
351  }
352 
353 
354  // remove outliers
355 
356  if (parameters.stdev > 0.0) {
357 
358  double chi2_old = evaluator(result, data.begin(), data.end());
359 
360  for ( ; ; ) {
361 
362  vector<hit_type>::iterator out = data.end();
363 
364  double xmax = 0.0;
365 
366  for (vector<hit_type>::iterator hit = data.begin(); hit != data.end(); ++hit) {
367 
368  const double x = fabs(hit->getValue() - estimator.getToA(result, *hit)) / hit->sigma;
369 
370  if (x > xmax) {
371  xmax = x;
372  out = hit;
373  }
374  }
375 
376  if (xmax > parameters.stdev) {
377 
378  DEBUG("Remove outlier " << out->getEKey() << ' ' << out->getLocation() << ' ' << FIXED(12,3) << xmax << "?" << flush);
379 
380  vector<hit_type>::iterator __end = data.end();
381 
382  iter_swap(out, --__end);
383 
384  result = estimator(data.begin(), __end);
385 
386  double chi2_new = evaluator(result, data.begin(), __end);
387 
388  if (chi2_old - chi2_new > parameters.stdev * parameters.stdev) {
389 
390  DEBUG(" yes " << FIXED(12,3) << (chi2_old - chi2_new) << endl);
391 
392  data.pop_back();
393 
394  chi2_old = chi2_new;
395 
396  continue;
397 
398  } else {
399 
400  DEBUG(" no " << endl);
401 
402  result = estimator(data.begin(), ++__end);
403  }
404  }
405 
406  break;
407  }
408 
409  chi2 = evaluator(result, data.begin(), data.end());
410  ndf = getWeight(data.begin(), data.end()) - result.getN();
411 
412  DEBUG("prefit:" << ' '
413  << FIXED(12,3) << chi2 << '/' << FIXED(6,1) << ndf << endl
414  << result << endl);
415 
416  for (vector<hit_type>::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
417  DEBUG("hit: " << *hit << ' ' << FIXED(9,3) << evaluator(result, *hit) << endl);
418  }
419  }
420 
421  for (vector<hit_type>::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
422  HB[hit->getID()]->Fill(geometry.getIndex(hit->getString()), hit->getFloor(), 1.0);
423  }
424 
425 
426  chi2 = numeric_limits<double>::max();
427 
428  switch (fit) {
429 
430  case linear_t:
431 
432  result = estimator(data.begin(), data.end());
433  chi2 = evaluator(result, data.begin(), data.end());
434 
435  hn.Fill(log10(1));
436  break;
437 
438  case simplex_t:
439 
440  simplex.value = result; // start value
441 
442  chi2 = simplex(data.begin(), data.end()) / simplex.estimator->getRho(1.0);
443  result = simplex.value;
444 
445  hn.Fill(log10(simplex.numberOfIterations));
446  break;
447 
448  case gandalf_t:
449 
450  gandalf.value = result; // start value
451 
452  chi2 = gandalf(data.begin(), data.end()) / gandalf.estimator->getRho(1.0);
453  result = gandalf.value;
454 
455  hn.Fill(log10(gandalf.numberOfIterations));
456  break;
457 
458  default:
459  break;
460  }
461 
462  ndf = getWeight(data.begin(), data.end()) - result.getN();
463 
464  DEBUG("result:" << ' '
465  << FIXED(12,3) << chi2 << '/' << FIXED(6,1) << ndf << endl
466  << result << endl);
467 
468  for (vector<hit_type>::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
469  DEBUG("hit: " << *hit << ' ' << FIXED(9,3) << evaluator(result, *hit) << endl);
470  }
471 
472  h0.Fill(chi2 / ndf);
473 
474  // store results
475 
476  if (chi2 / ndf <= parameters.chi2perNDF) {
477 
478  const JEvt evt = getEvt(JHead(oid,
479  data. begin()->getValue(),
480  data.rbegin()->getValue(),
481  data.size(),
482  result.getN(),
483  ndf,
484  chi2),
485  result);
486 
487  outputFile.put(evt);
488 
489  for (buffer_type::iterator i = p; i != q; ++i) {
490 
491  JEvent out(*i); // event with fitted times of emission
492 
493  const double toe = result.emitter[JEKey(i->getID(), distance(p,i))].t1;
494 
495  for (JEvent::iterator hit = out.begin(); hit != out.end(); ++hit) {
496  hit->setToE(toe);
497  }
498 
499  outputFile.put(out);
500  }
501 
502  counter[0] += 1;
503 
504  } else {
505 
506  WARNING(endl << "Event not written: " << FIXED(12,3) << chi2 << '/' << FIXED(6,1) << ndf << endl);
507 
508  counter[1] += 1;
509  }
510  }
511  }
512  }
513  STATUS(endl);
514 
515  STATUS("Number of events written / rejected: " << counter[0] << " / " << counter[1] << endl);
516 
517  JMultipleFileScanner<JMeta> io(inputFile);
518 
519  io >> outputFile;
520 
521  outputFile.put(h0);
522  outputFile.put(h1);
523  outputFile.put(hn);
524 
525  for (JManager<int, TH2D>::const_iterator i = HA.begin(); i != HA.end(); ++i) { outputFile.put(*(i->second)); }
526  for (JManager<int, TH2D>::const_iterator i = HB.begin(); i != HB.end(); ++i) { outputFile.put(*(i->second)); }
527 
528  outputFile.close();
529 }
Auxiliary class for ROOT I/O of application specific meta data.
Definition: JMeta.hh:70
const double xmax
Definition: JQuadrature.cc:24
Object writing to file.
Utility class to parse command line options.
Definition: JParser.hh:1517
size_t getN() const
Get number of fit parameters.
General exception.
Definition: JException.hh:23
#define WARNING(A)
Definition: JMessage.hh:65
double getValue(const JScale_t scale)
Get numerical value corresponding to scale.
Definition: JScale.hh:47
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
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
Template specialisation of fit function of acoustic model based on linear approximation.
Definition: JKatoomba.hh:249
Template specialisation of fit function of acoustic model based on JSimplex minimiser.
Definition: JKatoomba.hh:364
JEvt getEvt(const JHead &header, const JModel &model)
Get event.
*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
#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
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:446
V(JDAQEvent-JTriggerReprocessor)*1.0/(JDAQEvent+1.0e-10)
string outputFile
floor_range getRangeOfFloors(const JDetector &detector)
Get range of floors.
Acoustics hit.
size_t getNumberOfEmitters(T __begin, T __end)
Get number of emitters.
Linear fit.
Definition: JKatoomba.hh:47
Model for fit to acoustics data.
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.
Template specialisation of fit function of acoustic model based on JGandalf minimiser.
Definition: JKatoomba.hh:457
Acoustic transmission.
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
Acoustic event fit.
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
set_variable E_E log10(E_{fit}/E_{#mu})"
const std::string & getOID() const
Get detector identifier.
return result
Definition: JPolint.hh:764
JPosition3D getPosition(const Vec &pos)
Get position.
double getQ(const double D_m, const double f_kHz, const double d_m)
Get relative quality for given frequency at given distance.
Monte Carlo run header.
Definition: JHead.hh:1167
static const JStringCounter getNumberOfStrings
Function object to count unique strings.
Implementation for depth dependend velocity of sound.
#define FATAL(A)
Definition: JMessage.hh:67
then fatal Wrong number of arguments fi set_variable DETECTOR $argv[1] set_variable TRIPOD $argv[2] set_variable TX $argv[3] set_variable TY $argv[4] if[[!-f $DETECTOR]]
Definition: JFootprint.sh:28
z range($ZMAX-$ZMIN)< $MINIMAL_DZ." fi fi typeset -Z 4 STRING typeset -Z 2 FLOOR JPlot1D -f $
JACOUSTICS::JModel::emitter_type emitter
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.
Template specialisation of fit function of acoustic model based on JAbstractMinimiser minimiser...
Definition: JKatoomba.hh:188
then if[[!-f $DETECTOR]] then JDetector sh $DETECTOR fi cat $WORKDIR trigger_parameters txt<< EOFtrigger3DMuon.enabled=1;trigger3DMuon.numberOfHits=5;trigger3DMuon.gridAngle_deg=1;ctMin=0.0;TMaxLocal_ns=15.0;EOF set_variable TRIGGEREFFICIENCY_TRIGGERED_EVENTS_ONLY INPUT_FILES=() for((i=1;$i<=$NUMBER_OF_RUNS;++i));do JSirene.sh $DETECTOR $JPP_DATA/genhen.km3net_wpd_V2_0.evt.gz $WORKDIR/sirene_ ${i}.root JTriggerEfficiency.sh $DETECTOR $DETECTOR $WORKDIR/sirene_ ${i}.root $WORKDIR/trigger_efficiency_ ${i}.root $WORKDIR/trigger_parameters.txt $JPP_DATA/PMT_parameters.txt INPUT_FILES+=($WORKDIR/trigger_efficiency_ ${i}.root) done for ANGLE_DEG in $ANGLES_DEG[*];do set_variable SIGMA_NS 3.0 set_variable OUTLIERS 3 set_variable OUTPUT_FILE $WORKDIR/matrix\[${ANGLE_DEG}\deg\].root $JPP_DIR/examples/JReconstruction-f"$INPUT_FILES[*]"-o $OUTPUT_FILE-S ${SIGMA_NS}-A ${ANGLE_DEG}-O ${OUTLIERS}-d ${DEBUG}--!fiif[[$OPTION=="plot"]];then if((0));then for H1 in h0 h1;do JPlot1D-f"$WORKDIR/matrix["${^ANGLES_DEG}" deg].root:${H1}"-y"1 2e3"-Y-L TR-T""-\^"number of events [a.u.]"-> o chi2
Definition: JMatrixNZ.sh:106
Emitter key.
Definition: JEKey.hh:32
Acoustic event.
Custom probability density function of time-of-arrival.
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:73
JMEstimator * getMEstimator(const int type)
Get M-Estimator.
Definition: JMEstimator.hh:203
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.
Simplex fit.
Definition: JKatoomba.hh:48
Template definition of fit function of acoustic model.
Definition: JKatoomba.hh:180
int debug
debug level
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62
Gandalf fit.
Definition: JKatoomba.hh:49