Jpp  16.0.0
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
examples/JAcoustics/JKatoomba.cc File Reference

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

#include <string>
#include <iostream>
#include <iomanip>
#include <vector>
#include <set>
#include <map>
#include <algorithm>
#include "TRandom3.h"
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TH2D.h"
#include "JLang/JPredicate.hh"
#include "JLang/JSTDToolkit.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JTripod.hh"
#include "JDetector/JModule.hh"
#include "JDetector/JHydrophone.hh"
#include "JROOT/JManager.hh"
#include "JROOT/JRootToolkit.hh"
#include "JAcoustics/JSoundVelocity.hh"
#include "JAcoustics/JEmitter.hh"
#include "JAcoustics/JAcousticsToolkit.hh"
#include "JAcoustics/JAcousticsSupportkit.hh"
#include "JAcoustics/JHit.hh"
#include "JAcoustics/JFitParameters.hh"
#include "JAcoustics/JKatoomba.hh"
#include "JAcoustics/JTransmission_t.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 simulated acoustic data.

Author
mdejong

Definition in file examples/JAcoustics/JKatoomba.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 59 of file examples/JAcoustics/JKatoomba.cc.

60 {
61  using namespace std;
62  using namespace JPP;
63 
64  typedef JContainer< vector<JTripod> > tripods_container;
65  typedef JContainer< vector<JHydrophone> > hydrophones_container;
66  typedef JContainer< map<int, JWaveform> > waveform_container;
67  typedef JContainer< set<JTransmission_t> > disable_container;
68 
69  string detectorFile;
70  string outputFile;
71  int numberOfEvents;
72  map<int, int> numberOfPings;
73  JSoundVelocity V = getSoundVelocity; // default sound velocity
74  tripods_container tripods; // tripods
75  hydrophones_container hydrophones; // hydrophones
76  JTimeRange T_s(-1.0e-4, +1.0e-4); // time range of background [s]
77  JFitParameters parameters; // fit parameters
78  int fit; // type of fit
79  bool unify; // unify weighing of pings
80  disable_container disable; // disable tansmissions
81  waveform_container waveforms; // emitter power and frequency
82  UInt_t seed;
83  int debug;
84 
85  numberOfPings[WILDCARD] = NUMBER_OF_PINGS;
86  waveforms [WILDCARD] = EMITTER_WAVEFORM;
87 
88  try {
89 
90  JParser<> zap("Example application to test fit of model to simulated acoustic data.");
91 
92  zap['a'] = make_field(detectorFile);
93  zap['o'] = make_field(outputFile);
94  zap['n'] = make_field(numberOfEvents);
96  zap['N'] = make_field(numberOfPings, "number of pings per emitter") = JPARSER::initialised();
97  zap['V'] = make_field(V, "sound velocity") = JPARSER::initialised();
98  zap['T'] = make_field(tripods, "tripod data");
99  zap['H'] = make_field(hydrophones, "hydrophone data") = JPARSER::initialised();
100  zap['F'] = make_field(fit, "fit type") = linear_t, simplex_t, gandalf_t;
101  zap['u'] = make_field(unify, "unify weighing of pings");
102  zap['!'] = make_field(disable, "disable transmission") = JPARSER::initialised();
103  zap['W'] = make_field(waveforms, "emitter power and frequency") = JPARSER::initialised();
104  zap['S'] = make_field(seed) = 0;
105  zap['d'] = make_field(debug) = 1;
106 
107  zap(argc, argv);
108  }
109  catch(const exception &error) {
110  FATAL(error.what() << endl);
111  }
112 
113 
114  cout.tie(&cerr);
115 
116  gRandom->SetSeed(seed);
117 
118  if (numberOfEvents <= 0) { FATAL("Invalid number of events " << numberOfEvents << endl); }
119 
120 
122 
123  try {
124  load(detectorFile, detector);
125  }
126  catch(const JException& error) {
127  FATAL(error);
128  }
129 
130  vector<JEmitter> emitters;
131 
132  for (tripods_container::const_iterator i = tripods.begin(); i != tripods.end(); ++i) {
133  emitters.push_back(JEmitter(i->getID(), i->getUTMPosition() - detector.getUTMPosition()));
134  }
135 
136  if (detector.empty()) { FATAL("No modules in detector." << endl); }
137  if (emitters.empty()) { FATAL("No emitters in system." << endl); }
138 
139 
140  V.set(detector.getUTMZ()); // sound velocity at detector depth
141 
142  const double D_m = -detector.getUTMZ(); // depth [m]
143 
144 
145  const JGeometry geometry(detector, hydrophones);
146 
147  DEBUG(geometry);
148 
149 
150  JKatoomba<JEstimator> estimator(geometry, V);
151  JKatoomba<JAbstractMinimiser> evaluator(geometry, V);
152  JKatoomba<JSimplex> simplex (geometry, V);
153  JKatoomba<JGandalf> gandalf (geometry, V);
154 
155  simplex.estimator.reset(getMEstimator(parameters.mestimator));
156  gandalf.estimator.reset(getMEstimator(parameters.mestimator));
157 
158  simplex.debug = debug;
159  gandalf.debug = debug;
160 
162 
163 
164  typedef JHit<JPDFGauss> hit_type;
165 
166 
167  TH1D h0("cpu", NULL, 100, 1.0, 7.0);
168  TH1D h1("chi2/NDF", NULL, 100, 0.0, 5.0);
169 
170  JManager<int, TH1D> H1(new TH1D("emitter[%]", NULL, 100, -5.0e-5, +5.0e-5));
171  JManager<int, TH2D> H2(new TH2D("string[%]", NULL, 100, -1.0, +1.0, 100, -1.0, +1.0));
172 
173 
174  JTimer timer;
175 
176  for (int number_of_events = 0, count = 0; number_of_events != numberOfEvents; ++number_of_events) {
177 
178  STATUS("event: " << setw(10) << number_of_events << '\r'); DEBUG(endl);
179 
180  JModel model; // randomised model data
181 
182  for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) {
183 
184  if (model.string[module->getString()] == JMODEL::JString()) {
185 
186  model.string[module->getString()] = JMODEL::JString(gRandom->Uniform(-2.0e-2, +2.0e-2),
187  gRandom->Uniform(-2.0e-2, +2.0e-2));
188  }
189  }
190 
191  DEBUG("Model" << endl << model << endl);
192 
193  // set module positions according actual model data
194 
195  for (JDetector::iterator module = detector.begin(); module != detector.end(); ++module) {
196  if (geometry.hasLocation(module->getLocation())) {
197  module->set(estimator.detector[module->getString()].getPosition(model.string[module->getString()], module->getFloor()));
198  }
199  }
200 
201 
202  // generate hits
203 
204  int minimum_number_of_pings = numeric_limits<int>::max();
205 
206  for (vector<JEmitter>::const_iterator emitter = emitters.begin(); emitter != emitters.end(); ++emitter) {
207  minimum_number_of_pings = min(minimum_number_of_pings, getValue(numberOfPings, emitter->getID(), numberOfPings[WILDCARD]));
208  }
209 
210  vector<hit_type> data;
211 
212  for (vector<JEmitter>::const_iterator emitter = emitters.begin(); emitter != emitters.end(); ++emitter) {
213 
214  const int number_of_pings = getValue(numberOfPings, emitter->getID(), numberOfPings[WILDCARD]);
215 
216  const double signal = (unify ? (double) minimum_number_of_pings / (double) number_of_pings : 1.0);
217 
218  for (int ping_counter = 0; ping_counter != number_of_pings; ++ping_counter) {
219 
220  ++count;
221 
222  const double toe_s = ping_counter * 10.0 + gRandom->Uniform(-1.0, +1.0);
223 
224  model.emitter[emitter->getID()][count] = toe_s;
225 
226  for (JDetector::const_iterator module = detector.begin(); module != detector.end(); ++module) {
227 
228  if (disable.count(JTransmission_t(emitter->getID(), module->getID())) == 0) {
229 
230  if (geometry.hasLocation(module->getLocation())) {
231 
232  const JWaveform& waveform = getValue(waveforms, emitter->getID(), waveforms[WILDCARD]);
233 
234  const double d_m = getDistance(module->getPosition(), emitter->getPosition());
235  const double toa_s = toe_s + V.getTime(d_m, emitter->getZ(), module->getZ());
236  const double Q = waveform.getQ(D_m, d_m);
237 
238  if (Q >= parameters.Qmin) {
239 
240  double t1_s = toa_s;
241 
242  if (gRandom->Rndm() >= parameters.background)
243  t1_s = gRandom->Gaus(toa_s, parameters.sigma_s);
244  else
245  t1_s = gRandom->Uniform(toa_s + T_s.getLowerLimit(),
246  toa_s + T_s.getUpperLimit());
247 
248  const hit_type hit(*emitter,
249  count,
250  module->getLocation(),
251  JPDFGauss(t1_s, parameters.sigma_s, signal, parameters.background));
252 
253  DEBUG("hit: " << hit << ' ' << FIXED(7,1) << Q << endl);
254 
255  data.push_back(hit);
256  }
257  }
258  }
259  }
260  }
261  }
262 
263 
264  timer.reset();
265  timer.start();
266 
267  JModel result = estimator(data.begin(), data.end()); // pre-fit
268  double chi2 = numeric_limits<double>::max();
269 
270  switch (fit) {
271 
272  case linear_t:
273 
274  chi2 = evaluator(result, data.begin(), data.end());
275  break;
276 
277  case simplex_t:
278 
279  simplex.value = result; // start value
280 
281  chi2 = simplex(data.begin(), data.end()) / simplex.estimator->getRho(1.0);
282  result = simplex.value;
283  break;
284 
285  case gandalf_t:
286 
287  gandalf.value = result; // start value
288 
289  chi2 = gandalf(data.begin(), data.end()) / gandalf.estimator->getRho(1.0);
290  result = gandalf.value;
291  break;
292 
293  default:
294  break;
295  }
296 
297  timer.stop();
298 
299  double W = 0.0;
300 
301  for (vector<hit_type>::const_iterator hit = data.begin(); hit != data.end(); ++hit) {
302  W += hit->getWeight();
303  }
304 
305  const int ndf = data.size() - result.getN();
306 
307  DEBUG("Final values" << endl
308  << FIXED(9,3) << chi2 << '/' << ndf << endl
309  << result << endl);
310 
311 
312  h0.Fill(log10((double) timer.usec_wall));
313  h1.Fill(chi2 / (double) (W - result.getN()));
314 
315  for (JHashMap<int, JMODEL::JString>::const_iterator i = model.string.begin(); i != model.string.end(); ++i) {
316  H2[i->first]->Fill((i->second.tx - result.string [i->first].tx) * 1.0e3, // mrad
317  (i->second.ty - result.string [i->first].ty) * 1.0e3); // mrad
318  }
319 
320  for (JHashMap<JEKey, JMODEL::JEmitter>::const_iterator i = model.emitter.begin(); i != model.emitter.end(); ++i) {
321  H1[i->first.getID()]->Fill(i->second.t1 - result.emitter[i->first].t1);
322  }
323  }
324  STATUS(endl);
325 
326 
327  TFile out(outputFile.c_str(), "recreate");
328 
329  out << h0
330  << h1
331  << H1
332  << H2;
333 
334  out.Write();
335  out.Close();
336 }
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
Q(UTCMax_s-UTCMin_s)-livetime_s
double getValue(const JScale_t scale)
Get numerical value corresponding to scale.
Definition: JScale.hh:47
#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
*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
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
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:446
V(JDAQEvent-JTriggerReprocessor)*1.0/(JDAQEvent+1.0e-10)
string outputFile
Acoustics hit.
double getDistance(const JFirst_t &first, const JSecond_t &second)
Get distance between objects.
Linear fit.
Definition: JKatoomba.hh:45
Model for fit to acoustics data.
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
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:224
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
set_variable E_E log10(E_{fit}/E_{#mu})"
return result
Definition: JPolint.hh:743
Auxiliary class for CPU timing and usage.
Definition: JTimer.hh:32
double getQ(const double D_m, const double d_m) const
Get quality at given distance.
int debug
debug level
Definition: JSirene.cc:63
Implementation for depth dependend velocity of sound.
#define FATAL(A)
Definition: JMessage.hh:67
JACOUSTICS::JModel::string_type string
JACOUSTICS::JModel::emitter_type emitter
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
std::vector< int > count
Definition: JAlgorithm.hh:180
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
Custom probability density function of time-of-arrival.
container_type::const_iterator const_iterator
Definition: JHashMap.hh:86
JMEstimator * getMEstimator(const int type)
Get M-Estimator.
Definition: JMEstimator.hh:166
do set_variable DETECTOR_TXT $WORKDIR detector
then fatal Not enough tripods
Acoustic transmission identifier.
Utility class for emitter power and frequency.
Simplex fit.
Definition: JKatoomba.hh:46
Gandalf fit.
Definition: JKatoomba.hh:47