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 "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 "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 53 of file software/JAcoustics/JKatoomba.cc.

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