Jpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
JKatoomba.cc File Reference

Example application to test fit of model to acoustic data. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <vector>
#include <map>
#include <algorithm>
#include "TRandom3.h"
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TH2D.h"
#include "JLang/JPredicate.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JBase.hh"
#include "JDetector/JTripod.hh"
#include "JDetector/JModule.hh"
#include "JTools/JQuantile.hh"
#include "JGizmo/JManager.hh"
#include "JROOT/JRootToolkit.hh"
#include "JAcoustics/JSoundVelocity.hh"
#include "JAcoustics/JEmitter.hh"
#include "JAcoustics/JAcousticsToolkit.hh"
#include "JAcoustics/JHit.hh"
#include "JAcoustics/JKatoomba.hh"
#include "Jeep/JContainer.hh"
#include "Jeep/JTimer.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

Example application to test fit of model to acoustic data.

Author
mdejong

Definition in file JKatoomba.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 83 of file JKatoomba.cc.

84 {
85  using namespace std;
86  using namespace JPP;
87 
88  string detectorFile;
89  string outputFile;
90  int numberOfEvents;
91  int numberOfCycles;
92  JSoundVelocity V = getSoundVelocity; // default sound velocity
93  JContainer<JTripod> tripods; // tripods
94  JContainer<JBase> bases; // base modules
95  int mestimator;
96  UInt_t seed;
97  int debug;
98 
99  try {
100 
101  JParser<> zap("Example application to test fit of model to acoustic data.");
102 
103  zap['a'] = make_field(detectorFile);
104  zap['o'] = make_field(outputFile);
105  zap['n'] = make_field(numberOfEvents);
106  zap['N'] = make_field(numberOfCycles, "number of cycles per emitter") = 10;
107  zap['V'] = make_field(V, "sound velocity") = JPARSER::initialised();
108  zap['T'] = make_field(tripods, "tripod data");
109  zap['B'] = make_field(bases, "base module data") = JPARSER::initialised();
110  zap['M'] = make_field(mestimator) = EM_LINEAR, EM_LORENTZIAN, EM_NORMAL;
111  zap['S'] = make_field(seed) = 0;
112  zap['d'] = make_field(debug) = 1;
113 
114  zap(argc, argv);
115  }
116  catch(const exception &error) {
117  FATAL(error.what() << endl);
118  }
119 
120 
121  cout.tie(&cerr);
122 
123  gRandom->SetSeed(seed);
124 
125  const double sigma_s = 1.0e-5; // resolution of time-of_arrival [s]
126  const double P_bg = 1.0e-4; // probability of random background
127  const JTimeRange T_s(-1.0e-4, +1.0e-4); // time range of random background [s]
128 
129  if (numberOfEvents <= 0) { FATAL("Invalid number of events " << numberOfEvents << endl); }
130  if (numberOfCycles <= 0) { FATAL("Invalid number of cycles " << numberOfCycles << endl); }
131 
133 
134  try {
135  load(detectorFile, detector);
136  }
137  catch(const JException& error) {
138  FATAL(error);
139  }
140 
141  vector<JEmitter> emitters;
142 
143  for (JContainer<JTripod>::const_iterator i = tripods.begin(); i != tripods.end(); ++i) {
144  emitters.push_back(JEmitter(i->getID(), i->getUTMPosition() - detector.getUTMPosition()));
145  }
146 
147  if (detector.empty()) { FATAL("No modules in detector." << endl); }
148  if (emitters.empty()) { FATAL("No emitters in system." << endl); }
149 
150 
152 
153  {
154  map<int, vector<module_type> > buffer; // string -> module
155 
156  for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) {
157  buffer[module->getString()].push_back(module_type(module->getLocation(), module->getPosition()));
158  }
159 
160  for (map<int, vector<module_type> >::iterator i = buffer.begin(); i != buffer.end(); ++i) {
161 
162  JPosition3D p0; // reference position of string
163 
164  JContainer<JBase>::const_iterator p = find_if(bases.begin(), bases.end(), make_predicate(&JBase::getString, i->first));
165 
166  if (p != bases.end()) {
167 
168  p0 = p->getPosition();
169 
170  } else if (i->second.size() >= 2) {
171 
172  sort(i->second.begin(), i->second.end());
173 
174  const double z0 = TBARZ_M; // position of T-bar
175  const module_type& p1 = i->second[0];
176 
177  JVector3D slope;
178 
179  for (vector<module_type>::const_iterator p2 = i->second.begin(); ++p2 != i->second.end(); ) {
180  slope += (*p2 - p1) / p2->getDistance(p1);
181  }
182 
183  slope /= (double) (i->second.size() - 1);
184 
185  // linear extrapolation
186 
187  p0 = p1 + (z0 - p1.getZ()) * slope;
188 
189  } else {
190 
191  FATAL("Invalid data for string " << i->first << endl);
192 
193  continue;
194  }
195 
196  NOTICE("Position string "
197  << setw(4) << i->first << ' '
198  << FIXED(7,2) << p0.getX() << ' '
199  << FIXED(7,2) << p0.getY() << ' '
200  << FIXED(7,2) << p0.getZ() << ' '
201  << (p != bases.end() ? "user": "calculated") << endl);
202 
203  abc[i->first] = JGeometry::JString(p0, i->second.begin(), i->second.end());
204 
205  abc[i->first].hydrophone = JPosition3D(); // ?
206  }
207 
208  DEBUG(abc);
209  }
210 
211 
212  JKatoomba<JEstimator> estimator(abc, V);
213  JKatoomba<JAbstractMinimiser> evaluator(abc, V);
214 
215  evaluator.estimator.reset(getMEstimator(mestimator));
216 
217 
218  typedef JHit<JPDFGauss> hit_type;
219 
220 
221  TH1D h0("cpu", NULL, 100, 1.0, 7.0);
222  TH1D h1("chi2/NDF", NULL, 100, 0.0, 5.0);
223 
224  JManager<int, TH1D> H1(new TH1D("emitter[%]", NULL, 100, -2.0e-5, +2.0e-5));
225  JManager<int, TH2D> H2(new TH2D("string[%]", NULL, 100, -2.0e-4, +2.0e-4, 100, -2.0e-4, +2.0e-4));
226 
227 
228  JTimer timer;
229 
230  for (int number_of_events = 0, count = 0; number_of_events != numberOfEvents; ++number_of_events) {
231 
232  STATUS("event: " << setw(10) << number_of_events << '\r'); DEBUG(endl);
233 
234  JModel model; // randomised model data
235 
236  for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) {
237 
238  if (model.string[module->getString()] == JModel::JString()) {
239 
240  model.string[module->getString()] = JModel::JString(gRandom->Uniform(-1.0e-2, +1.0e-2),
241  gRandom->Uniform(-1.0e-2, +1.0e-2));
242  }
243  }
244 
245 
246  // set module positions according actual model data
247 
248  for (JDetector::iterator module = detector.begin(); module != detector.end(); ++module) {
249  module->set(estimator.detector[module->getString()].getPosition(model.string[module->getString()], module->getFloor()));
250  }
251 
252 
253  // generate hits
254 
255  vector<hit_type> data;
256 
257  for (vector<JEmitter>::const_iterator emitter = emitters.begin(); emitter != emitters.end(); ++emitter) {
258 
259  for (int number_of_cycles = 0; number_of_cycles != numberOfCycles; ++number_of_cycles) {
260 
261  ++count;
262 
263  const double toe_s = number_of_cycles * 10.0 + gRandom->Uniform(-1.0, +1.0);
264 
265  model.emitter[emitter->getID()][count] = toe_s;
266 
267  for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) {
268 
269  const double toa_s = toe_s + V.getTime(module->getDistance(emitter->getPosition()), emitter->getZ(), module->getZ());
270 
271  double t1_s = toa_s;
272 
273  if (gRandom->Rndm() >= P_bg)
274  t1_s = gRandom->Gaus(toa_s, sigma_s);
275  else
276  t1_s = gRandom->Uniform(toa_s + T_s.getLowerLimit(),
277  toa_s + T_s.getUpperLimit());
278 
279  data.push_back(hit_type(*emitter,
280  count,
281  module->getLocation(),
282  JPDFGauss(t1_s, sigma_s, P_bg)));
283  }
284  }
285  }
286 
287  DEBUG("Model" << endl << model << endl);
288 
289  for (vector<hit_type>::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
290  DEBUG("hit: " << *hit << endl);
291  }
292 
293 
294  timer.reset();
295  timer.start();
296 
297  JModel result = estimator(data.begin(), data.end()); // pre-fit
298  double chi2 = evaluator(result, data.begin(), data.end());
299 
300  timer.stop();
301 
302 
303  DEBUG("Final values" << endl
304  << FIXED(9,3) << chi2 << '/'
305  << data.size() - model.getN() << endl
306  << result << endl);
307 
308 
309  h0.Fill(log10((double) timer.usec_wall));
310  h1.Fill(2.0 * chi2 / (double) (data.size() - model.getN()));
311 
312  for (JHashMap<int, JModel::JString>::const_iterator i = model.string.begin(); i != model.string.end(); ++i) {
313  H2[i->first]->Fill(i->second.tx - result.string [i->first].tx,
314  i->second.ty - result.string [i->first].ty);
315  }
316 
317  for (JHashMap<JEKey, JModel::JEmitter>::const_iterator i = model.emitter.begin(); i != model.emitter.end(); ++i) {
318  H1[i->first.getID()]->Fill(i->second.t1 - result.emitter[i->first].t1);
319  }
320  }
321  STATUS(endl);
322 
323 
324  TFile out(outputFile.c_str(), "recreate");
325 
326  out << h0
327  << h1
328  << H1
329  << H2;
330 
331  out.Write();
332  out.Close();
333 }
Utility class to parse command line options.
Definition: JParser.hh:1493
size_t getN() const
Get number of fit parameters.
General exception.
Definition: JException.hh:23
JPredicate< JResult_t T::*, JComparison::eq > make_predicate(JResult_t T::*member, const JResult_t value)
Helper method to create predicate for data member.
Definition: JPredicate.hh:128
TPaveText * p1
#define STATUS(A)
Definition: JMessage.hh:63
Detector data structure.
Definition: JDetector.hh:80
std::vector< T >::const_iterator const_iterator
Definition: JContainer.hh:32
Template specialisation of fit function of acoustic model based on linear approximation.
Definition: JKatoomba.hh:153
JPosition3D hydrophone
Hydrophone.
Definition: JGeometry.hh:304
then for HISTOGRAM in h0 h1
Definition: JMatrixNZ.sh:69
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition: JParser.hh:63
Auxiliary data structure for floating point format specification.
Definition: JPrint.hh:481
string outputFile
double getDistance(const JVector3D &pos) const
Get distance to point.
Definition: JVector3D.hh:269
Acoustics hit.
JACOUSTICS::JModel::emitter_map emitter
Auxiliary class to manage set of compatible ROOT objects (e.g.
Definition: JManager.hh:42
Model for fit to acoustics data.
JACOUSTICS::JModel::string_map string
static const JSoundVelocity getSoundVelocity(1541.0,-17.0e-3,-2000.0)
Function object for velocity of sound.
Detector file.
Definition: JHead.hh:130
Data structure for vector in three dimensions.
Definition: JVector3D.hh:33
Acoustic emitter.
Definition: JEmitter.hh:27
Auxiliary container for detector elements.
Definition: JContainer.hh:29
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1954
return result
Definition: JPolint.hh:695
Auxiliary class for CPU timing and usage.
Definition: JTimer.hh:32
#define NOTICE(A)
Definition: JMessage.hh:64
double getY() const
Get y position.
Definition: JVector3D.hh:103
void load(const JString &file_name, JDetector &detector)
Load detector from input file.
int debug
debug level
Definition: JSirene.cc:61
Implementation for velocity of sound.
#define FATAL(A)
Definition: JMessage.hh:67
std::vector< int > count
Definition: JAlgorithm.hh:184
Template specialisation of fit function of acoustic model based on JAbstractMinimiser minimiser...
Definition: JKatoomba.hh:94
double getX() const
Get x position.
Definition: JVector3D.hh:93
Custom probability density function of time-of-arrival.
Data structure for position in three dimensions.
Definition: JPosition3D.hh:35
container_type::const_iterator const_iterator
Definition: JHashMap.hh:85
JMEstimator * getMEstimator(const int type)
Get M-Estimator.
Definition: JMEstimator.hh:166
static const double TBARZ_M
T-bar position relative to seabed [m].
double getZ() const
Get z position.
Definition: JVector3D.hh:114
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62