Jpp  15.0.0
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 <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/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 55 of file software/JAcoustics/JKatoomba.cc.

56 {
57  using namespace std;
58  using namespace JPP;
59 
60  typedef JContainer< vector<JTripod> > tripods_container;
61  typedef JContainer< vector<JHydrophone> > hydrophones_container;
62  typedef JContainer< set<JTransmission_t> > disable_container;
63 
64  JMultipleFileScanner<> inputFile;
66  string detectorFile;
67  JLimit_t& numberOfEvents = inputFile.getLimit();
68  JSoundVelocity V = getSoundVelocity; // default sound velocity
69  tripods_container tripods; // tripods
70  hydrophones_container hydrophones; // hydrophones
71  JFitParameters parameters; // fit parameters
72  int fit; // type of fit
73  bool unify; // unify weighing of pings
74  disable_container disable; // disable tansmissions
75  int debug;
76 
77  try {
78 
79  JParser<> zap("Application to fit position calibration model to acoustic data.");
80 
81  zap['f'] = make_field(inputFile, "output of JAcousticEventBuilder[.sh]");
82  zap['o'] = make_field(outputFile);
83  zap['n'] = make_field(numberOfEvents) = JLimit::max();
84  zap['a'] = make_field(detectorFile);
86  zap['V'] = make_field(V, "sound velocity") = JPARSER::initialised();
87  zap['T'] = make_field(tripods, "tripod data");
88  zap['H'] = make_field(hydrophones, "hydrophone data") = JPARSER::initialised();
89  zap['M'] = make_field(getMechanics, "mechanics data") = JPARSER::initialised();
90  zap['F'] = make_field(fit, "fit type") = linear_t, simplex_t, gandalf_t;
91  zap['u'] = make_field(unify, "unify weighing of pings");
92  zap['!'] = make_field(disable, "disable transmission") = JPARSER::initialised();
93  zap['d'] = make_field(debug) = 1;
94 
95  zap(argc, argv);
96  }
97  catch(const exception &error) {
98  FATAL(error.what() << endl);
99  }
100 
101 
102  cout.tie(&cerr);
103 
105 
106  try {
107  load(detectorFile, detector);
108  }
109  catch(const JException& error) {
110  FATAL(error);
111  }
112 
114 
115  JHashMap<int, JLocation> receivers;
116  JHashMap<int, JEmitter> emitters;
117 
118  for (JDetector::const_iterator i = detector.begin(); i != detector.end(); ++i) {
119  receivers[i->getID()] = i->getLocation();
120  }
121 
122  for (tripods_container::const_iterator i = tripods.begin(); i != tripods.end(); ++i) {
123  emitters[i->getID()] = JEmitter(i->getID(),
124  i->getUTMPosition() - detector.getUTMPosition());
125  }
126 
127  V.set(detector.getUTMZ()); // sound velocity at detector depth
128 
129  JGeometry geometry(detector, hydrophones);
130 
131  DEBUG(geometry);
132 
133 
134  JKatoomba<JEstimator> estimator(geometry, V);
135  JKatoomba<JAbstractMinimiser> evaluator(geometry, V);
136  JKatoomba<JSimplex> simplex (geometry, V);
137  JKatoomba<JGandalf> gandalf (geometry, V);
138 
139  evaluator.estimator.reset(getMEstimator(parameters.mestimator));
140  simplex .estimator.reset(getMEstimator(parameters.mestimator));
141  gandalf .estimator.reset(getMEstimator(parameters.mestimator));
142 
143  simplex.debug = debug;
144  gandalf.debug = debug;
145 
146  if (parameters.fixStrings) {
147  JKatoomba_t::setOption(false);
148  }
149 
151 
152  typedef JHit<JPDFGauss> hit_type;
153 
154  TH1D h0("chi2/NDF", NULL, 50000, 0.0, 1000.0);
155  TH1D h1("h1", NULL, 51, -0.5, 50.5);
156 
157  JManager<int, TH2D> HA(new TH2D("ha[%]", NULL,
159  range.getLength(), range.getLowerLimit() - 0.5, range.getLowerLimit() + 0.5));
160 
161  JManager<int, TH2D> HB(new TH2D("hb[%]", NULL,
163  range.getLength(), range.getLowerLimit() - 0.5, range.getLowerLimit() + 0.5));
164 
165  for (Int_t i = 1; i <= HA->GetXaxis()->GetNbins(); ++i) {
166  HA->GetXaxis()->SetBinLabel(i, MAKE_CSTRING(geometry.at(i-1).first));
167  HB->GetXaxis()->SetBinLabel(i, MAKE_CSTRING(geometry.at(i-1).first));
168  }
169 
170  // input data
171 
172  string oid; // detector identifier
173 
174  vector<JEvent> zbuf;
175 
176  JLimit_t limit = inputFile.getLimit();
177  counter_type number_of_events = 0;
178 
179  for (JMultipleFileScanner<>::const_iterator file = inputFile.begin(); file != inputFile.end(); ++file) {
180 
182 
183  limit.setLowerLimit(limit.getLowerLimit() - in.skip(limit.getLowerLimit()));
184 
185  const size_t N = zbuf.size();
186 
187  for ( ; in.hasNext() && number_of_events != limit; ++number_of_events) {
188 
189  STATUS("input " << setw(8) << number_of_events << '\r'); DEBUG(endl);
190 
191  const JEvent* evt = in.next();
192 
193  if (oid == "") {
194  oid = evt->getOID();
195  }
196 
197  if (oid != evt->getOID()) { // consistency check
198  FATAL("Invalid detector identifier " << evt->getOID() << " != " << oid << endl);
199  }
200 
201  zbuf.push_back(*evt);
202  }
203 
204  vector<JEvent>::iterator q = zbuf.begin();
205 
206  advance(q, N);
207 
208  sort(q, zbuf.end()); // sort according first time-of-emission
209 
210  if (q != zbuf.begin() && q != zbuf.end()) {
211 
213 
214  --p;
215 
216  if (p->empty() || q->empty() || *q < *p) {
217  inplace_merge(zbuf.begin(), q, zbuf.end());
218  }
219  }
220  }
221  STATUS(endl);
222 
223 
224  outputFile.open();
225 
226  outputFile.put(JMeta(argc, argv));
227 
228  int counter[] = { 0, 0 };
229 
230  for (vector<JEvent>::const_iterator p = zbuf.begin(), q; p != zbuf.end(); p = q) {
231 
232  STATUS("output " << setw(8) << distance(zbuf.cbegin(), p) << '\r'); DEBUG(endl);
233 
234  for (q = p; ++q != zbuf.end() && q->begin()->getToE() <= p->rbegin()->getToE() + parameters.Tmax_s; ) {}
235 
236  if (getNumberOfEmitters(p,q) >= parameters.Nmin) {
237 
238  h1.Fill((double) getNumberOfEmitters(p,q));
239 
240  map<int, int> numberOfPings;
241 
242  for (vector<JEvent>::const_iterator i = p; i != q; ++i) {
243  numberOfPings[i->getID()] += 1;
244  }
245 
246  for (map<int, int>::const_iterator i = numberOfPings.begin(); i != numberOfPings.end(); ++i) {
247  DEBUG("Number of pings " << setw(2) << i->first << ' ' << setw(3) << i->second << endl);
248  }
249 
250  int minimum_number_of_pings = numeric_limits<int>::max();
251 
252  for (map<int, int>::const_iterator i = numberOfPings.begin(); i != numberOfPings.end(); ++i) {
253  minimum_number_of_pings = min(minimum_number_of_pings, i->second);
254  }
255 
256  vector<hit_type> data;
257 
258  for (vector<JEvent>::const_iterator i = p; i != q; ++i) {
259 
260  const JEmitter& emitter = emitters[i->getID()];
261 
262  const double signal = (unify ? (double) minimum_number_of_pings / (double) numberOfPings[i->getID()] : 1.0);
263 
264  // select transmission with highest quality
265 
266  map<int, vector<JTransmission> > buffer;
267 
268  for (JEvent::const_iterator hit = i->begin(); hit != i->end(); ++hit) {
269  if (disable.count(JTransmission_t(i->getID(), hit->getID())) == 0 &&
270  disable.count(JTransmission_t(-1, hit->getID())) == 0) {
271  buffer[hit->getID()].push_back(*hit);
272  }
273  }
274 
275  for (map<int, vector<JTransmission> >::iterator ps = buffer.begin(); ps != buffer.end(); ++ps) {
276 
277  if (receivers.has(ps->first) && geometry.hasLocation(receivers[ps->first])) {
278 
279  sort(ps->second.begin(), ps->second.end(), make_comparator(&JTransmission::getQ, JComparison::gt()));
280 
281  if (ps->second.begin()->getQ() >= parameters.Qmin) {
282 
283  data.push_back(hit_type(emitter,
284  i->getCounter(),
285  receivers[ps->first],
286  JPDFGauss(ps->second.begin()->getToA(), parameters.sigma_s, signal, parameters.background)));
287  }
288  }
289  }
290  }
291 
292  for (vector<hit_type>::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
293  DEBUG("hit: " << *hit << endl);
294  }
295 
296 
297 
298  JModel result = estimator(data.begin(), data.end()); // pre-fit
299 
300  DEBUG("prefit:" << endl
301  << result << endl);
302 
303 
304  for (vector<hit_type>::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
305  HA[hit->getID()]->Fill(geometry.getIndex(hit->getString()), hit->getFloor(), 1.0);
306  }
307 
308 
309  // remove outliers
310 
311  if (parameters.stdev > 0.0) {
312 
313  double chi2_old = evaluator(result, data.begin(), data.end());
314 
315  for ( ; ; ) {
316 
317  vector<hit_type>::iterator out = data.end();
318 
319  double xmax = 0.0;
320 
321  for (vector<hit_type>::iterator hit = data.begin(); hit != data.end(); ++hit) {
322 
323  const double x = fabs(hit->getValue() - estimator.getToA(result, *hit)) / hit->sigma;
324 
325  if (x > xmax) {
326  xmax = x;
327  out = hit;
328  }
329  }
330 
331  if (xmax > parameters.stdev) {
332 
333  vector<hit_type>::iterator __end = data.end();
334 
335  iter_swap(out, --__end);
336 
337  result = estimator(data.begin(), __end);
338 
339  double chi2_new = evaluator(result, data.begin(), __end);
340 
341  if (chi2_old - chi2_new > parameters.stdev * parameters.stdev) {
342 
343  DEBUG("Remove outlier " << __end->getLocation() << ' ' << xmax << endl);
344 
345  data.pop_back();
346 
347  chi2_old = chi2_new;
348 
349  continue;
350  }
351  }
352 
353  break;
354  }
355  }
356 
357  for (vector<hit_type>::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
358  DEBUG("hit: " << *hit << endl);
359  }
360 
361  DEBUG("prefit:" << endl
362  << result << endl);
363 
364  for (vector<hit_type>::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
365  HB[hit->getID()]->Fill(geometry.getIndex(hit->getString()), hit->getFloor(), 1.0);
366  }
367 
368 
369  double chi2 = numeric_limits<double>::max();
370 
371  switch (fit) {
372 
373  case linear_t:
374 
375  result = estimator(data.begin(), data.end());
376  chi2 = evaluator(result, data.begin(), data.end());
377  break;
378 
379  case simplex_t:
380 
381  simplex.value = result; // start value
382 
383  chi2 = simplex(data.begin(), data.end()) / simplex.estimator->getRho(1.0);
384  result = simplex.value;
385  break;
386 
387  case gandalf_t:
388 
389  gandalf.value = result; // start value
390 
391  chi2 = gandalf(data.begin(), data.end()) / gandalf.estimator->getRho(1.0);
392  result = gandalf.value;
393  break;
394 
395  default:
396  break;
397  }
398 
399  double W = 0.0;
400 
401  for (vector<hit_type>::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
402  W += hit->getWeight();
403  }
404 
405  const int ndf = data.size() - result.getN();
406 
407  DEBUG("result:" << endl
408  << FIXED(9,3) << chi2 << '/' << ndf << endl
409  << result << endl);
410 
411  h0.Fill(chi2 / (W - result.getN()));
412 
413  // store results
414 
415  if (chi2 / ndf <= parameters.chi2perNDF) {
416 
417  const JEvt evt = getEvt(JHead(oid, data.begin()->getValue(), data.rbegin()->getValue(), ndf, (W - result.getN()), chi2), result);
418 
419  outputFile.put(evt);
420 
421  for (vector<JEvent>::const_iterator i = p; i != q; ++i) {
422 
423  JEvent out(*i); // event with fitted times of emission
424 
425  const double toe = result.emitter[JEKey(i->getID(), i->getCounter())].t1;
426 
427  for (JEvent::iterator hit = out.begin(); hit != out.end(); ++hit) {
428  hit->setToE(toe);
429  }
430 
431  outputFile.put(out);
432  }
433 
434  counter[0] += 1;
435 
436  } else {
437 
438  WARNING(endl << "Event not written: " << chi2 << '/' << ndf << endl);
439 
440  counter[1] += 1;
441  }
442  }
443  }
444  STATUS(endl);
445 
446  STATUS("Number of events written / rejected: " << counter[0] << " / " << counter[1] << endl);
447 
448  JMultipleFileScanner<JMeta> io(inputFile);
449 
450  io >> outputFile;
451 
452  outputFile.put(h0);
453  outputFile.put(h1);
454 
455  for (JManager<int, TH2D>::const_iterator i = HA.begin(); i != HA.end(); ++i) { outputFile.put(*(i->second)); }
456  for (JManager<int, TH2D>::const_iterator i = HB.begin(); i != HB.end(); ++i) { outputFile.put(*(i->second)); }
457 
458  outputFile.close();
459 }
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:1500
size_t getN() const
Get number of fit parameters.
General exception.
Definition: JException.hh:23
#define WARNING(A)
Definition: JMessage.hh:65
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:285
Template specialisation of fit function of acoustic model based on JSimplex minimiser.
Definition: JKatoomba.hh:385
then JShowerPostfit f $INPUT_FILE o $OUTPUT_FILE N
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:151
then for HISTOGRAM in h0 h1
Definition: JMatrixNZ.sh:71
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition: JParser.hh:66
Long64_t counter_type
Type definition for counter.
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:446
string outputFile
floor_range getRangeOfFloors(const JDetector &detector)
Get range of floors.
void inplace_merge(T __begin, const size_t N, const size_t *delimiter)
Merge multiple sorted ranges.
Definition: JMergeSort.cc:29
Acoustics hit.
size_t getNumberOfEmitters(T __begin, T __end)
Get number of emitters.
Linear fit.
Definition: JKatoomba.hh:45
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:466
Acoustic transmission.
do for((RUN=${RANGE%%-*};$RUN<=${RANGE##*-};RUN+=1))
Auxiliary class to manage set of compatible ROOT objects (e.g. histograms) using unique keys...
Definition: JManager.hh:43
Detector file.
Definition: JHead.hh:196
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:1961
const std::string & getOID() const
Get detector identifier.
return result
Definition: JPolint.hh:727
int debug
debug level
Definition: JSirene.cc:63
double getQ(const double D_m, const double f_kHz, const double d_m)
Get relative quality for given frequency at given distance.
counter_type advance(counter_type &counter, const counter_type value, const counter_type limit=std::numeric_limits< counter_type >::max())
Advance counter.
Monte Carlo run header.
Definition: JHead.hh:1113
static const JStringCounter getNumberOfStrings
Function object to count unique strings.
Implementation for velocity of sound.
#define FATAL(A)
Definition: JMessage.hh:67
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.
Simple fit method based on Powell&#39;s algorithm, see reference: Numerical Recipes in C++...
Definition: JSimplex.hh:42
Template specialisation of fit function of acoustic model based on JAbstractMinimiser minimiser...
Definition: JKatoomba.hh:226
then usage $script< input_file >< detector_file > fi set_variable OUTPUT_DIR set_variable SELECTOR JDAQTimesliceL1 set_variable DEBUG case set_variable DEBUG
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:166
do set_variable DETECTOR_TXT $WORKDIR detector
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 source JAcoustics sh $DETECTOR_ID CHECK_EXIT_CODE typeset A TRIPODS get_tripods $WORKDIR tripod txt TRIPODS for EMITTER in
Definition: JCanberra.sh:41
Acoustic transmission identifier.
Simplex fit.
Definition: JKatoomba.hh:46
do if[[!-f $ACOUSTICS_WORKDIR/${KEY}.txt]]
Definition: JAcoustics.sh:39
then usage $script< string identifier >< detectorfile > event file(toashort file)+" "\nNote that the event files and toashort files should be one-to-one related." fi if (( $
Gandalf fit.
Definition: JKatoomba.hh:47