Jpp - 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 fit position calibration model to acoustic data. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <vector>
#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 "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 "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 fit position calibration model to acoustic data.

Author
mdejong

Definition in file software/JAcoustics/JKatoomba.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

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

52 {
53  using namespace std;
54  using namespace JPP;
55 
56  typedef JContainer< vector<JTripod> > tripods_container;
57  typedef JContainer< vector<JHydrophone> > hydrophones_container;
58 
61  string detectorFile;
62  JLimit_t& numberOfEvents = inputFile.getLimit();
63  JSoundVelocity V = getSoundVelocity; // default sound velocity
64  tripods_container tripods; // tripods
65  hydrophones_container hydrophones; // hydrophones
66  JFitParameters parameters; // fit parameters
67  int fit; // type of fit
68  bool unique; // one ping per cycle
69  int debug;
70 
71  try {
72 
73  JParser<> zap("Application to fit position calibration model to acoustic data.");
74 
75  zap['f'] = make_field(inputFile, "output of JAcousticEventBuilder[.sh]");
76  zap['o'] = make_field(outputFile);
77  zap['n'] = make_field(numberOfEvents) = JLimit::max();
78  zap['a'] = make_field(detectorFile);
80  zap['V'] = make_field(V, "sound velocity") = JPARSER::initialised();
81  zap['T'] = make_field(tripods, "tripod data");
82  zap['H'] = make_field(hydrophones, "hydrophone data") = JPARSER::initialised();
83  zap['M'] = make_field(getMechanics, "mechanics data") = JPARSER::initialised();
84  zap['F'] = make_field(fit, "fit type") = linear_t, simplex_t, gandalf_t;
85  zap['u'] = make_field(unique, "one ping per cycle");
86  zap['d'] = make_field(debug) = 1;
87 
88  zap(argc, argv);
89  }
90  catch(const exception &error) {
91  FATAL(error.what() << endl);
92  }
93 
94 
95  cout.tie(&cerr);
96 
98 
99  try {
100  load(detectorFile, detector);
101  }
102  catch(const JException& error) {
103  FATAL(error);
104  }
105 
106 
107  JHashMap<int, JLocation> receivers;
108  JHashMap<int, JEmitter> emitters;
109 
110  for (JDetector::const_iterator i = detector.begin(); i != detector.end(); ++i) {
111  receivers[i->getID()] = i->getLocation();
112  }
113 
114  for (tripods_container::const_iterator i = tripods.begin(); i != tripods.end(); ++i) {
115  emitters[i->getID()] = JEmitter(i->getID(),
116  i->getUTMPosition() - detector.getUTMPosition());
117  }
118 
119 
120  JGeometry geometry(detector, hydrophones);
121 
122  DEBUG(geometry);
123 
124 
125  JKatoomba<JEstimator> estimator(geometry, V[detector.getUTMZ()]);
126  JKatoomba<JAbstractMinimiser> evaluator(geometry, V[detector.getUTMZ()]);
127  JKatoomba<JSimplex> simplex (geometry, V[detector.getUTMZ()]);
128  JKatoomba<JGandalf> gandalf (geometry, V[detector.getUTMZ()]);
129 
130  evaluator.estimator.reset(getMEstimator(parameters.mestimator));
131  simplex .estimator.reset(getMEstimator(parameters.mestimator));
132  gandalf .estimator.reset(getMEstimator(parameters.mestimator));
133 
134  simplex.debug = debug;
135  gandalf.debug = debug;
136 
138 
139  typedef JHit<JPDFGauss> hit_type;
140 
141  TH1D h0("chi2/NDF", NULL, 5000, 0.0, 100.0);
142 
143  TH2D ha("ha", NULL,
146 
147  TH2D hb("hb", NULL,
150 
151  for (Int_t i = 1; i <= ha.GetXaxis()->GetNbins(); ++i) {
152  ha.GetXaxis()->SetBinLabel(i, MAKE_CSTRING(geometry.at(i-1).first));
153  hb.GetXaxis()->SetBinLabel(i, MAKE_CSTRING(geometry.at(i-1).first));
154  }
155 
156  // input data
157 
158  string oid; // detector identifier
159 
160  vector<JEvent> zbuf;
161 
162  while (inputFile.hasNext()) {
163 
164  STATUS("input " << setw(6) << inputFile.getCounter() << '\r'); DEBUG(endl);
165 
166  const JEvent* evt = inputFile.next();
167 
168  if (oid == "") {
169  oid = evt->getOID();
170  }
171 
172  if (oid != evt->getOID()) { // consistency check
173  FATAL("Invalid detector identifier " << evt->getOID() << " != " << oid << endl);
174  }
175 
176  zbuf.push_back(*evt);
177  }
178  STATUS(endl);
179 
180  sort(zbuf.begin(), zbuf.end());
181 
182 
183  outputFile.open();
184 
185  outputFile.put(JMeta(argc, argv));
186 
187 
188  for (vector<JEvent>::const_iterator p = zbuf.begin(), q; p != zbuf.end(); p = q) {
189 
190  STATUS("output " << setw(6) << distance(zbuf.cbegin(), p) << '\r'); DEBUG(endl);
191 
192  for (q = p; ++q != zbuf.end() && q->begin()->getToE() <= p->rbegin()->getToE() + parameters.Tmax_s; ) {}
193 
194  if (getNumberOfEmitters(p,q) >= parameters.Nmin) {
195 
196  map<int, int> numberOfPings;
197 
198  for (vector<JEvent>::const_iterator i = p; i != q; ++i) {
199  numberOfPings[i->getID()] += 1;
200  }
201 
202  int weight = numeric_limits<int>::max();
203 
204  for (map<int, int>::const_iterator i = numberOfPings.begin(); i != numberOfPings.end(); ++i) {
205  DEBUG("Number of pings " << setw(2) << i->first << ' ' << setw(3) << i->second << endl);
206  }
207 
208  for (map<int, int>::const_iterator i = numberOfPings.begin(); i != numberOfPings.end(); ++i) {
209  if (i->second < weight) {
210  weight = i->second;
211  }
212  }
213 
214  vector<hit_type> data;
215 
216  for (vector<JEvent>::const_iterator i = p; i != q; ++i) {
217 
218  {
219  const JEmitter& emitter = emitters[i->getID()];
220 
221  const double signal = (unique ? (double) weight / (double) numberOfPings[i->getID()] : 1.0);
222 
223  // select first transmission
224 
225  map<int, vector<JTransmission> > buffer;
226 
227  for (JEvent::const_iterator hit = i->begin(); hit != i->end(); ++hit) {
228  buffer[hit->getID()].push_back(*hit);
229  }
230 
231  for (map<int, vector<JTransmission> >::iterator ps = buffer.begin(); ps != buffer.end(); ++ps) {
232 
233  if (receivers.has(ps->first) && geometry.hasLocation(receivers[ps->first])) {
234 
235  sort(ps->second.begin(), ps->second.end(), make_comparator(&JTransmission::getToA));
236 
237  data.push_back(hit_type(emitter,
238  i->getCounter(),
239  receivers[ps->first],
240  JPDFGauss(ps->second.begin()->getToA(), parameters.sigma_s, signal, parameters.background)));
241  }
242  }
243  }
244  }
245 
246  for (vector<hit_type>::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
247  DEBUG("hit: " << *hit << endl);
248  }
249 
250 
251 
252  JModel result = estimator(data.begin(), data.end()); // pre-fit
253 
254  DEBUG("prefit:" << endl
255  << result << endl);
256 
257 
258  for (vector<hit_type>::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
259  ha.Fill(geometry.getIndex(hit->getString()), hit->getFloor(), 1.0);
260  }
261 
262 
263  // remove outliers
264 
265  if (parameters.stdev > 0.0) {
266 
267  double chi2_old = evaluator(result, data.begin(), data.end());
268 
269  for ( ; ; ) {
270 
271  vector<hit_type>::iterator out = data.end();
272 
273  double xmax = 0.0;
274 
275  for (vector<hit_type>::iterator hit = data.begin(); hit != data.end(); ++hit) {
276 
277  const double x = fabs(hit->getValue() - estimator.getToA(result, *hit)) / hit->sigma;
278 
279  if (x > xmax) {
280  xmax = x;
281  out = hit;
282  }
283  }
284 
285  if (xmax > parameters.stdev) {
286 
287  vector<hit_type>::iterator __end = data.end();
288 
289  iter_swap(out, --__end);
290 
291  result = estimator(data.begin(), __end);
292 
293  double chi2_new = evaluator(result, data.begin(), __end);
294 
295  if (chi2_old - chi2_new > parameters.stdev * parameters.stdev) {
296 
297  DEBUG("Remove outlier " << __end->getLocation() << ' ' << xmax << endl);
298 
299  data.pop_back();
300 
301  chi2_old = chi2_new;
302 
303  continue;
304  }
305  }
306 
307  break;
308  }
309  }
310 
311  for (vector<hit_type>::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
312  DEBUG("hit: " << *hit << endl);
313  }
314 
315  DEBUG("prefit:" << endl
316  << result << endl);
317 
318  for (vector<hit_type>::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
319  hb.Fill(geometry.getIndex(hit->getString()), hit->getFloor(), 1.0);
320  }
321 
322 
323  double chi2 = numeric_limits<double>::max();
324 
325  switch (fit) {
326 
327  case linear_t:
328 
329  result = estimator(data.begin(), data.end());
330  chi2 = evaluator(result, data.begin(), data.end());
331  break;
332 
333  case simplex_t:
334 
335  simplex.value = result; // start value
336 
337  chi2 = simplex(data.begin(), data.end()) / simplex.estimator->getRho(1.0);
338  result = simplex.value;
339  break;
340 
341  case gandalf_t:
342 
343  gandalf.value = result; // start value
344 
345  chi2 = gandalf(data.begin(), data.end()) / gandalf.estimator->getRho(1.0);
346  result = gandalf.value;
347  break;
348 
349  default:
350  break;
351  }
352 
353  double W = 0.0;
354 
355  for (vector<hit_type>::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
356  W += hit->getWeight();
357  }
358 
359  const int ndf = data.size() - result.getN();
360 
361  DEBUG("result:" << endl
362  << FIXED(9,3) << chi2 << '/' << ndf << endl
363  << result << endl);
364 
365  h0.Fill(chi2 / (W - result.getN()));
366 
367  // store results
368 
369  const JEvt evt = getEvt(JHead(oid, data.begin()->getValue(), data.rbegin()->getValue(), ndf, chi2), result);
370 
371  outputFile.put(evt);
372 
373  for (vector<JEvent>::const_iterator i = p; i != q; ++i) {
374 
375  JEvent out(*i);
376 
377  const double toe = result.emitter[JEKey(i->getID(), i->getCounter())].t1;
378 
379  for (JEvent::iterator hit = out.begin(); hit != out.end(); ++hit) {
380  hit->setToE(toe);
381  }
382 
383  outputFile.put(out);
384  }
385  }
386  }
387  STATUS(endl);
388 
389  JMultipleFileScanner<JMeta> io(inputFile);
390 
391  io >> outputFile;
392 
393  outputFile.put(h0);
394  outputFile.put(ha);
395  outputFile.put(hb);
396 
397  outputFile.close();
398 }
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
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:71
#define STATUS(A)
Definition: JMessage.hh:63
Detector data structure.
Definition: JDetector.hh:80
Template specialisation of fit function of acoustic model based on linear approximation.
Definition: JKatoomba.hh:250
Template specialisation of fit function of acoustic model based on JSimplex minimiser.
Definition: JKatoomba.hh:345
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
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition: JParser.hh:66
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:446
string outputFile
Acoustics hit.
size_t getNumberOfEmitters(T __begin, T __end)
Get number of emitters.
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:426
Acoustic transmission.
do for((RUN=${RANGE%%-*};$RUN<=${RANGE##*-};RUN+=1))
Detector file.
Definition: JHead.hh:196
Acoustic event fit.
Acoustic emitter.
Definition: JEmitter.hh:27
JLANG::JSharedPointer< JMEstimator > estimator
M-Estimator function.
Definition: JKatoomba.hh:94
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 getNumberOfFloors(const JDetector &detector)
Get number of floors.
int debug
debug level
Definition: JSirene.cc:63
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
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.
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:191
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
std::vector< double > weight
Definition: JAlgorithm.hh:428
do if[[!-f $ACOUSTICS_WORKDIR/${KEY}.txt]]
Definition: JAcoustics.sh:38
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62